Linear Algebra

This chapter describes functions for solving linear systems. The library provides linear algebra operations which operate directly on the gsl_vector and gsl_matrix objects. These routines use the standard algorithms from Golub & Van Loan’s Matrix Computations with Level-1 and Level-2 BLAS calls for efficiency.

The functions described in this chapter are declared in the header file gsl_linalg.h.

LU Decomposition

A general M-by-N matrix A has an LU decomposition

P A = L U

where P is an M-by-M permutation matrix, L is M-by-\min(M,N) and U is \min(M,N)-by-N. For square matrices, L is a lower unit triangular matrix and U is upper triangular. For M > N, L is a unit lower trapezoidal matrix of size M-by-N. For M < N, U is upper trapezoidal of size M-by-N. For square matrices this decomposition can be used to convert the linear system A x = b into a pair of triangular systems (L y = P b, U x = y), which can be solved by forward and back-substitution. Note that the LU decomposition is valid for singular matrices.

int gsl_linalg_LU_decomp(gsl_matrix *A, gsl_permutation *p, int *signum)
int gsl_linalg_complex_LU_decomp(gsl_matrix_complex *A, gsl_permutation *p, int *signum)

These functions factorize the matrix A into the LU decomposition PA = LU. On output the diagonal and upper triangular (or trapezoidal) part of the input matrix A contain the matrix U. The lower triangular (or trapezoidal) part of the input matrix (excluding the diagonal) contains L. The diagonal elements of L are unity, and are not stored.

The permutation matrix P is encoded in the permutation p on output. The j-th column of the matrix P is given by the k-th column of the identity matrix, where k = p_j the j-th element of the permutation vector. The sign of the permutation is given by signum. It has the value (-1)^n, where n is the number of interchanges in the permutation.

The algorithm used in the decomposition is Gaussian Elimination with partial pivoting (Golub & Van Loan, Matrix Computations, Algorithm 3.4.1), combined with a recursive algorithm based on Level 3 BLAS (Peise and Bientinesi, 2016).

The functions return GSL_SUCCESS for non-singular matrices. If the matrix is singular, the factorization is still completed, and the functions return an integer k \in [1, MIN(M,N)] such that U(k,k) = 0. In this case, U is singular and the decomposition should not be used to solve linear systems.

int gsl_linalg_LU_solve(const gsl_matrix *LU, const gsl_permutation *p, const gsl_vector *b, gsl_vector *x)
int gsl_linalg_complex_LU_solve(const gsl_matrix_complex *LU, const gsl_permutation *p, const gsl_vector_complex *b, gsl_vector_complex *x)

These functions solve the square system A x = b using the LU decomposition of A into (LU, p) given by gsl_linalg_LU_decomp() or gsl_linalg_complex_LU_decomp() as input.

int gsl_linalg_LU_svx(const gsl_matrix *LU, const gsl_permutation *p, gsl_vector *x)
int gsl_linalg_complex_LU_svx(const gsl_matrix_complex *LU, const gsl_permutation *p, gsl_vector_complex *x)

These functions solve the square system A x = b in-place using the precomputed LU decomposition of A into (LU, p). On input x should contain the right-hand side b, which is replaced by the solution on output.

int gsl_linalg_LU_refine(const gsl_matrix *A, const gsl_matrix *LU, const gsl_permutation *p, const gsl_vector *b, gsl_vector *x, gsl_vector *work)
int gsl_linalg_complex_LU_refine(const gsl_matrix_complex *A, const gsl_matrix_complex *LU, const gsl_permutation *p, const gsl_vector_complex *b, gsl_vector_complex *x, gsl_vector_complex *work)

These functions apply an iterative improvement to x, the solution of A x = b, from the precomputed LU decomposition of A into (LU, p). Additional workspace of length N is required in work.

int gsl_linalg_LU_invert(const gsl_matrix *LU, const gsl_permutation *p, gsl_matrix *inverse)
int gsl_linalg_complex_LU_invert(const gsl_matrix_complex *LU, const gsl_permutation *p, gsl_matrix_complex *inverse)

These functions compute the inverse of a matrix A from its LU decomposition (LU, p), storing the result in the matrix inverse. The inverse is computed by computing the inverses U^{-1}, L^{-1} and finally forming the product A^{-1} = U^{-1} L^{-1} P. Each step is based on Level 3 BLAS calls.

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).

int gsl_linalg_LU_invx(gsl_matrix *LU, const gsl_permutation *p)
int gsl_linalg_complex_LU_invx(gsl_matrix_complex *LU, const gsl_permutation *p)

These functions compute the inverse of a matrix A from its LU decomposition (LU, p), storing the result in-place in the matrix LU. The inverse is computed by computing the inverses U^{-1}, L^{-1} and finally forming the product A^{-1} = U^{-1} L^{-1} P. Each step is based on Level 3 BLAS calls.

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).

double gsl_linalg_LU_det(gsl_matrix *LU, int signum)
gsl_complex gsl_linalg_complex_LU_det(gsl_matrix_complex *LU, int signum)

These functions compute the determinant of a matrix A from its LU decomposition, LU. The determinant is computed as the product of the diagonal elements of U and the sign of the row permutation signum.

double gsl_linalg_LU_lndet(gsl_matrix *LU)
double gsl_linalg_complex_LU_lndet(gsl_matrix_complex *LU)

These functions compute the logarithm of the absolute value of the determinant of a matrix A, \ln|\det(A)|, from its LU decomposition, LU. This function may be useful if the direct computation of the determinant would overflow or underflow.

int gsl_linalg_LU_sgndet(gsl_matrix *LU, int signum)
gsl_complex gsl_linalg_complex_LU_sgndet(gsl_matrix_complex *LU, int signum)

These functions compute the sign or phase factor of the determinant of a matrix A, \det(A)/|\det(A)|, from its LU decomposition, LU.

QR Decomposition

A general rectangular M-by-N matrix A has a QR decomposition into the product of a unitary M-by-M square matrix Q (where Q^{\dagger} Q = I) and an M-by-N right-triangular matrix R,

A = Q R

This decomposition can be used to convert the square linear system A x = b into the triangular system R x = Q^{\dagger} b, which can be solved by back-substitution. Another use of the QR decomposition is to compute an orthonormal basis for a set of vectors. The first N columns of Q form an orthonormal basis for the range of A, ran(A), when A has full column rank.

When M > N, the bottom M - N rows of R are zero, and so A can be naturally partioned as

A = \begin{pmatrix} Q_1 & Q_2 \end{pmatrix} \begin{pmatrix} R_1 \\ 0 \end{pmatrix} = Q_1 R_1

where R_1 is N-by-N upper triangular, Q_1 is M-by-N, and Q_2 is M-by-(M-N). Q_1 R_1 is sometimes called the thin or reduced QR decomposition. The solution of the least squares problem \min_x || b - A x ||^2 when A has full rank is:

x = R_1^{-1} c_1

where c_1 is the first N elements of Q^{\dagger} b. If A is rank deficient, see QR Decomposition with Column Pivoting and Complete Orthogonal Decomposition.

GSL offers two interfaces for the QR decomposition. The first proceeds by zeroing out columns below the diagonal of A, one column at a time using Householder transforms. In this method, the factor Q is represented as a product of Householder reflectors:

Q = H_n \cdots H_2 H_1

where each H_i = I - \tau_i v_i v_i^{\dagger} for a scalar \tau_i and column vector v_i. In this method, functions which compute the full matrix Q or apply Q^{\dagger} to a right hand side vector operate by applying the Householder matrices one at a time using Level 2 BLAS.

The second interface is based on a Level 3 BLAS block recursive algorithm developed by Elmroth and Gustavson. In this case, Q is written in block form as

Q = I - V T V^{\dagger}

where V is an M-by-N matrix of the column vectors v_i and T is an N-by-N upper triangular matrix, whose diagonal elements are the \tau_i. Computing the full T, while requiring more flops than the Level 2 approach, offers the advantage that all standard operations can take advantage of cache efficient Level 3 BLAS operations, and so this method often performs faster than the Level 2 approach. The functions for the recursive block algorithm have a _r suffix, and it is recommended to use this interface for performance critical applications.

int gsl_linalg_QR_decomp_r(gsl_matrix *A, gsl_matrix *T)
int gsl_linalg_complex_QR_decomp_r(gsl_matrix_complex *A, gsl_matrix_complex *T)

These functions factor the M-by-N matrix A into the QR decomposition A = Q R using the recursive Level 3 BLAS algorithm of Elmroth and Gustavson. On output the diagonal and upper triangular part of A contain the matrix R. The N-by-N matrix T stores the upper triangular factor appearing in Q. The matrix Q is given by Q = I - V T V^{\dagger}, where the elements below the diagonal of A contain the columns of V on output.

This algorithm requires M \ge N and performs best for “tall-skinny” matrices, i.e. M \gg N.

int gsl_linalg_QR_solve_r(const gsl_matrix *QR, const gsl_matrix *T, const gsl_vector *b, gsl_vector *x)
int gsl_linalg_complex_QR_solve_r(const gsl_matrix_complex *QR, const gsl_matrix_complex *T, const gsl_vector_complex *b, gsl_vector_complex *x)

These functions solve the square system A x = b using the QR decomposition of A held in (QR, T). The least-squares solution for rectangular systems can be found using gsl_linalg_QR_lssolve_r() or gsl_linalg_complex_QR_lssolve_r().

int gsl_linalg_QR_lssolve_r(const gsl_matrix *QR, const gsl_matrix *T, const gsl_vector *b, gsl_vector *x, gsl_vector *work)
int gsl_linalg_complex_QR_lssolve_r(const gsl_matrix_complex *QR, const gsl_matrix_complex *T, const gsl_vector_complex *b, gsl_vector_complex *x, gsl_vector_complex *work)

These functions find the least squares solution to the overdetermined system A x = b, where the matrix A has more rows than columns. The least squares solution minimizes the Euclidean norm of the residual, ||b - Ax||. The routine requires as input the QR decomposition of A into (QR, T) given by gsl_linalg_QR_decomp_r() or gsl_linalg_complex_QR_decomp_r().

The parameter x is of length M. The solution x is returned in the first N rows of x, i.e. x = x[0], x[1], ..., x[N-1]. The last M - N rows of x contain a vector whose norm is equal to the residual norm || b - A x ||. This similar to the behavior of LAPACK DGELS. Additional workspace of length N is required in work.

int gsl_linalg_QR_lssolvem_r(const gsl_matrix *QR, const gsl_matrix *T, const gsl_matrix *B, gsl_matrix *X, gsl_matrix *work)
int gsl_linalg_complex_QR_lssolvem_r(const gsl_matrix_complex *QR, const gsl_matrix_complex *T, const gsl_matrix_complex *B, gsl_matrix_complex *X, gsl_matrix_complex *work)

These functions find the least squares solutions to the overdetermined systems A x_k = b_k where the matrix A has more rows than columns. The least squares solution minimizes the Euclidean norm of the residual, ||b_k - A x_k||. The routine requires as input the QR decomposition of A into (QR, T) given by gsl_linalg_QR_decomp_r() or gsl_linalg_complex_QR_decomp_r(). The right hand side b_k is provided in column k of the input B, while the solution x_k is stored in the first N rows of column k of the output X.

The parameters X and B are of size M-by-nrhs. The last M - N rows of X contain vectors whose norm is equal to the residual norm || b_k - A x_k ||. This similar to the behavior of LAPACK DGELS. Additional workspace of length N-by-nrhs is required in work.

int gsl_linalg_QR_QTvec_r(const gsl_matrix *QR, const gsl_matrix *T, gsl_vector *v, gsl_vector *work)
int gsl_linalg_complex_QR_QHvec_r(const gsl_matrix_complex *QR, const gsl_matrix_complex *T, gsl_vector_complex *v, gsl_vector_complex *work)

These functions apply the matrix Q^T (or Q^{\dagger}) encoded in the decomposition (QR, T) to the vector v, storing the result Q^T v (or Q^{\dagger} v) in v. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q. Additional workspace of size N is required in work.

int gsl_linalg_QR_QTmat_r(const gsl_matrix *QR, const gsl_matrix *T, gsl_matrix *B, gsl_matrix *work)
int gsl_linalg_complex_QR_QHmat_r(const gsl_matrix_complex *QR, const gsl_matrix_complex *T, gsl_matrix_complex *B, gsl_matrix_complex *work)

This function applies the matrix Q^T (or Q^{\dagger}) encoded in the decomposition (QR, T) to the M-by-K matrix B, storing the result Q^T B (or Q^{\dagger} B) in B. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q. Additional workspace of size N-by-K is required in work.

int gsl_linalg_QR_unpack_r(const gsl_matrix *QR, const gsl_matrix *T, gsl_matrix *Q, gsl_matrix *R)
int gsl_linalg_complex_QR_unpack_r(const gsl_matrix_complex *QR, const gsl_matrix_complex *T, gsl_matrix_complex *Q, gsl_matrix_complex *R)

These functions unpack the encoded QR decomposition (QR, T) as output from gsl_linalg_QR_decomp_r() or gsl_linalg_complex_QR_decomp_r() into the matrices Q and R, where Q is M-by-M and R is N-by-N. Note that the full R matrix is M-by-N, however the lower trapezoidal portion is zero, so only the upper triangular factor is stored.

int gsl_linalg_QR_rcond(const gsl_matrix *QR, double *rcond, gsl_vector *work)

This function estimates the reciprocal condition number (using the 1-norm) of the R factor, stored in the upper triangle of QR. The reciprocal condition number estimate, defined as 1 / (||R||_1 \cdot ||R^{-1}||_1), is stored in rcond. Additional workspace of size 3 N is required in work.

Level 2 Interface

The functions below are for the slower Level 2 interface to the QR decomposition. It is recommended to use these functions only for M < N, since the Level 3 interface above performs much faster for M \ge N.

int gsl_linalg_QR_decomp(gsl_matrix *A, gsl_vector *tau)
int gsl_linalg_complex_QR_decomp(gsl_matrix_complex *A, gsl_vector_complex *tau)

These functions factor the M-by-N matrix A into the QR decomposition A = Q R. On output the diagonal and upper triangular part of the input matrix contain the matrix R. The vector tau and the columns of the lower triangular part of the matrix A contain the Householder coefficients and Householder vectors which encode the orthogonal matrix Q. The vector tau must be of length N. The matrix Q is related to these components by the product of k=min(M,N) reflector matrices, Q = H_k ... H_2 H_1 where H_i = I - \tau_i v_i v_i^{\dagger} and v_i is the Householder vector v_i = (0,...,1,A(i+1,i),A(i+2,i),...,A(m,i)). This is the same storage scheme as used by LAPACK.

The algorithm used to perform the decomposition is Householder QR (Golub & Van Loan, “Matrix Computations”, Algorithm 5.2.1).

int gsl_linalg_QR_solve(const gsl_matrix *QR, const gsl_vector *tau, const gsl_vector *b, gsl_vector *x)
int gsl_linalg_complex_QR_solve(const gsl_matrix_complex *QR, const gsl_vector_complex *tau, const gsl_vector_complex *b, gsl_vector_complex *x)

These functions solve the square system A x = b using the QR decomposition of A held in (QR, tau). The least-squares solution for rectangular systems can be found using gsl_linalg_QR_lssolve().

int gsl_linalg_QR_svx(const gsl_matrix *QR, const gsl_vector *tau, gsl_vector *x)
int gsl_linalg_complex_QR_svx(const gsl_matrix_complex *QR, const gsl_vector_complex *tau, gsl_vector_complex *x)

These functions solve the square system A x = b in-place using the QR decomposition of A held in (QR, tau). On input x should contain the right-hand side b, which is replaced by the solution on output.

int gsl_linalg_QR_lssolve(const gsl_matrix *QR, const gsl_vector *tau, const gsl_vector *b, gsl_vector *x, gsl_vector *residual)
int gsl_linalg_complex_QR_lssolve(const gsl_matrix_complex *QR, const gsl_vector_complex *tau, const gsl_vector_complex *b, gsl_vector_complex *x, gsl_vector_complex *residual)

These functions find the least squares solution to the overdetermined system A x = b where the matrix A has more rows than columns. The least squares solution minimizes the Euclidean norm of the residual, ||Ax - b||.The routine requires as input the QR decomposition of A into (QR, tau) given by gsl_linalg_QR_decomp() or gsl_linalg_complex_QR_decomp(). The solution is returned in x. The residual is computed as a by-product and stored in residual.

int gsl_linalg_QR_QTvec(const gsl_matrix *QR, const gsl_vector *tau, gsl_vector *v)
int gsl_linalg_complex_QR_QHvec(const gsl_matrix_complex *QR, const gsl_vector_complex *tau, gsl_vector_complex *v)

These functions apply the matrix Q^T (or Q^{\dagger}) encoded in the decomposition (QR, tau) to the vector v, storing the result Q^T v (or Q^{\dagger} v) in v. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q^T (or Q^{\dagger}).

int gsl_linalg_QR_Qvec(const gsl_matrix *QR, const gsl_vector *tau, gsl_vector *v)
int gsl_linalg_complex_QR_Qvec(const gsl_matrix_complex *QR, const gsl_vector_complex *tau, gsl_vector_complex *v)

These functions apply the matrix Q encoded in the decomposition (QR, tau) to the vector v, storing the result Q v in v. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q.

int gsl_linalg_QR_QTmat(const gsl_matrix *QR, const gsl_vector *tau, gsl_matrix *B)

This function applies the matrix Q^T encoded in the decomposition (QR, tau) to the M-by-K matrix B, storing the result Q^T B in B. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q^T.

int gsl_linalg_QR_Rsolve(const gsl_matrix *QR, const gsl_vector *b, gsl_vector *x)

This function solves the triangular system R x = b for x. It may be useful if the product b' = Q^T b has already been computed using gsl_linalg_QR_QTvec().

int gsl_linalg_QR_Rsvx(const gsl_matrix *QR, gsl_vector *x)

This function solves the triangular system R x = b for x in-place. On input x should contain the right-hand side b and is replaced by the solution on output. This function may be useful if the product b' = Q^T b has already been computed using gsl_linalg_QR_QTvec().

int gsl_linalg_QR_unpack(const gsl_matrix *QR, const gsl_vector *tau, gsl_matrix *Q, gsl_matrix *R)

This function unpacks the encoded QR decomposition (QR, tau) into the matrices Q and R, where Q is M-by-M and R is M-by-N.

int gsl_linalg_QR_QRsolve(gsl_matrix *Q, gsl_matrix *R, const gsl_vector *b, gsl_vector *x)

This function solves the system R x = Q^T b for x. It can be used when the QR decomposition of a matrix is available in unpacked form as (Q, R).

int gsl_linalg_QR_update(gsl_matrix *Q, gsl_matrix *R, gsl_vector *w, const gsl_vector *v)

This function performs a rank-1 update w v^T of the QR decomposition (Q, R). The update is given by Q'R' = Q (R + w v^T) where the output matrices Q and R are also orthogonal and right triangular. Note that w is destroyed by the update.

int gsl_linalg_R_solve(const gsl_matrix *R, const gsl_vector *b, gsl_vector *x)

This function solves the triangular system R x = b for the N-by-N matrix R.

int gsl_linalg_R_svx(const gsl_matrix *R, gsl_vector *x)

This function solves the triangular system R x = b in-place. On input x should contain the right-hand side b, which is replaced by the solution on output.

Triangle on Top of Rectangle

This section provides routines for computing the QR decomposition of the specialized matrix

\begin{pmatrix} U \\ A \end{pmatrix} = Q R

where U is an N-by-N upper triangular matrix, and A is an M-by-N dense matrix. This type of matrix arises, for example, in the sequential TSQR algorithm. The Elmroth and Gustavson algorithm is used to efficiently factor this matrix. Due to the upper triangular factor, the Q matrix takes the form

Q = I - V T V^T

with

V = \begin{pmatrix} I \\ Y \end{pmatrix}

and Y is dense and of the same dimensions as A.

int gsl_linalg_QR_UR_decomp(gsl_matrix *U, gsl_matrix *A, gsl_matrix *T)

This function computes the QR decomposition of the matrix (U ; A), where U is N-by-N upper triangular and A is M-by-N dense. On output, U is replaced by the R factor, and A is replaced by Y. The N-by-N upper triangular block reflector is stored in T on output.

int gsl_linalg_QR_UR_lssolve(const gsl_matrix *R, const gsl_matrix *Y, const gsl_matrix *T, const gsl_vector *b, gsl_vector *x, gsl_vector *work)

This function finds the least squares solution to the overdetermined system,

\min_x \left| \left| b - \begin{pmatrix} U \\ A \end{pmatrix} x \right| \right|^2

where U is a N-by-N upper triangular matrix, and A is a M-by-N dense matrix. The routine requires as input the QR decomposition of (U; A) into (R, Y, T) given by gsl_linalg_QR_UR_decomp(). The parameter x is of length N+M. The solution x is returned in the first N rows of x, i.e. x = x[0], x[1], ..., x[N-1]. The last M rows of x contain a vector whose norm is equal to the residual norm || b - (U; A) x ||. This similar to the behavior of LAPACK DGELS. Additional workspace of length N is required in work.

int gsl_linalg_QR_UR_lssvx(const gsl_matrix *R, const gsl_matrix *Y, const gsl_matrix *T, gsl_vector *x, gsl_vector *work)

This function finds the least squares solution to the overdetermined system,

\min_x \left| \left| b - \begin{pmatrix} U \\ A \end{pmatrix} x \right| \right|^2

in-place, where U is a N-by-N upper triangular matrix, and A is a M-by-N dense matrix. The routine requires as input the QR decomposition of (U; A) into (R, Y, T) given by gsl_linalg_QR_UR_decomp(). The parameter x is of length N+M and contains the right hand side vector b on input. The solution x is returned in the first N rows of x, i.e. x = x[0], x[1], ..., x[N-1]. The last M rows of x contain a vector whose norm is equal to the residual norm || b - (U; A) x ||. This similar to the behavior of LAPACK DGELS. Additional workspace of length N is required in work.

int gsl_linalg_QR_UR_QTvec(const gsl_matrix *Y, const gsl_matrix *T, gsl_vector *b, gsl_vector *work)

This function computes Q^T b using the decomposition (Y, T) previously computed by gsl_linalg_QR_UR_decomp(). On input, b contains the length N+M vector b, and on output it will contain Q^T b. Additional workspace of length N is required in work.

Triangle on Top of Triangle

This section provides routines for computing the QR decomposition of the specialized matrix

\begin{pmatrix} U_1 \\ U_2 \end{pmatrix} = Q R

where U_1,U_2 are N-by-N upper triangular matrices. The Elmroth and Gustavson algorithm is used to efficiently factor this matrix. The Q matrix takes the form

Q = I - V T V^T

with

V = \begin{pmatrix} I \\ Y \end{pmatrix}

and Y is N-by-N upper triangular.

int gsl_linalg_QR_UU_decomp(gsl_matrix *U1, gsl_matrix *U2, gsl_matrix *T)

This function computes the QR decomposition of the matrix (U_1 ; U_2), where U_1,U_2 are N-by-N upper triangular. On output, U1 is replaced by the R factor, and U2 is replaced by Y. The N-by-N upper triangular block reflector is stored in T on output.

int gsl_linalg_QR_UU_lssolve(const gsl_matrix *R, const gsl_matrix *Y, const gsl_matrix *T, const gsl_vector *b, gsl_vector *x, gsl_vector *work)

This function finds the least squares solution to the overdetermined system,

\min_x \left| \left| b - \begin{pmatrix} U_1 \\ U_2 \end{pmatrix} x \right| \right|^2

where U_1,U_2 are N-by-N upper triangular matrices. The routine requires as input the QR decomposition of (U_1; U_2) into (R, Y, T) given by gsl_linalg_QR_UU_decomp(). The parameter x is of length 2N. The solution x is returned in the first N rows of x, i.e. x = x[0], x[1], ..., x[N-1]. The last N rows of x contain a vector whose norm is equal to the residual norm || b - (U_1; U_2) x ||. This similar to the behavior of LAPACK DGELS. Additional workspace of length N is required in work.

int gsl_linalg_QR_UU_lssvx(const gsl_matrix *R, const gsl_matrix *Y, const gsl_matrix *T, gsl_vector *x, gsl_vector *work)

This function finds the least squares solution to the overdetermined system,

\min_x \left| \left| b - \begin{pmatrix} U_1 \\ U_2 \end{pmatrix} x \right| \right|^2

in-place, where U_1,U_2 are N-by-N upper triangular matrices. The routine requires as input the QR decomposition of (U_1; U_2) into (R, Y, T) given by gsl_linalg_QR_UU_decomp(). The parameter x is of length 2N and contains the right hand side vector b on input. The solution x is returned in the first N rows of x, i.e. x = x[0], x[1], ..., x[N-1]. The last N rows of x contain a vector whose norm is equal to the residual norm || b - (U_1; U_2) x ||. This similar to the behavior of LAPACK DGELS. Additional workspace of length N is required in work.

int gsl_linalg_QR_UU_QTvec(const gsl_matrix *Y, const gsl_matrix *T, gsl_vector *b, gsl_vector *work)

This function computes Q^T b using the decomposition (Y, T) previously computed by gsl_linalg_QR_UU_decomp(). On input, b contains the vector b, and on output it will contain Q^T b. Additional workspace of length N is required in work.

Triangle on Top of Trapezoidal

This section provides routines for computing the QR decomposition of the specialized matrix

\begin{pmatrix} U \\ A \end{pmatrix} = Q R

where U is an N-by-N upper triangular matrix, and A is an M-by-N upper trapezoidal matrix with M \ge N. A has the structure,

A = \begin{pmatrix} A_d \\ A_u \end{pmatrix}

where A_d is (M-N)-by-N dense, and A_u is N-by-N upper triangular. The Elmroth and Gustavson algorithm is used to efficiently factor this matrix. The Q matrix takes the form

Q = I - V T V^T

with

V = \begin{pmatrix} I \\ Y \end{pmatrix}

and Y is upper trapezoidal and of the same dimensions as A.

int gsl_linalg_QR_UZ_decomp(gsl_matrix *U, gsl_matrix *A, gsl_matrix *T)

This function computes the QR decomposition of the matrix (U ; A), where U is N-by-N upper triangular and A is M-by-N upper trapezoidal. On output, U is replaced by the R factor, and A is replaced by Y. The N-by-N upper triangular block reflector is stored in T on output.

Triangle on Top of Diagonal

This section provides routines for computing the QR decomposition of the specialized matrix

\begin{pmatrix} U \\ D \end{pmatrix} = Q R

where U is an N-by-N upper triangular matrix and D is an N-by-N diagonal matrix. This type of matrix arises in regularized least squares problems. The Elmroth and Gustavson algorithm is used to efficiently factor this matrix. The Q matrix takes the form

Q = I - V T V^T

with

V = \begin{pmatrix} I \\ Y \end{pmatrix}

and Y is N-by-N upper triangular.

int gsl_linalg_QR_UD_decomp(gsl_matrix *U, const gsl_vector *D, gsl_matrix *Y, gsl_matrix *T)

This function computes the QR decomposition of the matrix (U ; D), where U is N-by-N upper triangular and D is N-by-N diagonal. On output, U is replaced by the R factor and Y is stored in Y. The N-by-N upper triangular block reflector is stored in T on output.

int gsl_linalg_QR_UD_lssolve(const gsl_matrix *R, const gsl_matrix *Y, const gsl_matrix *T, const gsl_vector *b, gsl_vector *x, gsl_vector *work)

This function finds the least squares solution to the overdetermined system,

\min_x \left| \left| b - \begin{pmatrix} U \\ D \end{pmatrix} x \right| \right|^2

where U is N-by-N upper triangular and D is N-by-N diagonal. The routine requires as input the QR decomposition of (U; D) into (R, Y, T) given by gsl_linalg_QR_UD_decomp(). The parameter x is of length 2N. The solution x is returned in the first N rows of x, i.e. x = x[0], x[1], ..., x[N-1]. The last N rows of x contain a vector whose norm is equal to the residual norm || b - (U; D) x ||. This similar to the behavior of LAPACK DGELS. Additional workspace of length N is required in work.

int gsl_linalg_QR_UD_lssvx(const gsl_matrix *R, const gsl_matrix *Y, const gsl_matrix *T, gsl_vector *x, gsl_vector *work)

This function finds the least squares solution to the overdetermined system,

\min_x \left| \left| b - \begin{pmatrix} U \\ D \end{pmatrix} x \right| \right|^2

in-place, where U is N-by-N upper triangular and D is N-by-N diagonal. The routine requires as input the QR decomposition of (U; D) into (R, Y, T) given by gsl_linalg_QR_UD_decomp(). The parameter x is of length 2N and contains the right hand side vector b on input. The solution x is returned in the first N rows of x, i.e. x = x[0], x[1], ..., x[N-1]. The last N rows of x contain a vector whose norm is equal to the residual norm || b - (U; D) x ||. This similar to the behavior of LAPACK DGELS. Additional workspace of length N is required in work.

int gsl_linalg_QR_UD_QTvec(const gsl_matrix *Y, const gsl_matrix *T, gsl_vector *b, gsl_vector *work)

This function computes Q^T b using the decomposition (Y, T) previously computed by gsl_linalg_QR_UD_decomp(). On input, b contains the vector b, and on output it will contain Q^T b. Additional workspace of length N is required in work.

QR Decomposition with Column Pivoting

The QR decomposition of an M-by-N matrix A can be extended to the rank deficient case by introducing a column permutation P,

A P = Q R

The first r columns of Q form an orthonormal basis for the range of A for a matrix with column rank r. This decomposition can also be used to convert the square linear system A x = b into the triangular system R y = Q^T b, x = P y, which can be solved by back-substitution and permutation. We denote the QR decomposition with column pivoting by QRP^T since A = Q R P^T. When A is rank deficient with r = {\rm rank}(A), the matrix R can be partitioned as

R = \left(
\begin{matrix}
  R_{11} & R_{12} \\
  0 & R_{22}
\end{matrix}
\right) \approx
\left(
\begin{matrix}
  R_{11} & R_{12} \\
  0 & 0
\end{matrix}
\right)

where R_{11} is r-by-r and nonsingular. In this case, a basic least squares solution for the overdetermined system A x = b can be obtained as

x = P \left(
\begin{matrix}
  R_{11}^{-1} c_1 \\
  0
\end{matrix}
\right)

where c_1 consists of the first r elements of Q^T b. This basic solution is not guaranteed to be the minimum norm solution unless R_{12} = 0 (see Complete Orthogonal Decomposition).

int gsl_linalg_QRPT_decomp(gsl_matrix *A, gsl_vector *tau, gsl_permutation *p, int *signum, gsl_vector *norm)

This function factorizes the M-by-N matrix A into the QRP^T decomposition A = Q R P^T. On output the diagonal and upper triangular part of the input matrix contain the matrix R. The permutation matrix P is stored in the permutation p. The sign of the permutation is given by signum. It has the value (-1)^n, where n is the number of interchanges in the permutation. The vector tau and the columns of the lower triangular part of the matrix A contain the Householder coefficients and vectors which encode the orthogonal matrix Q. The vector tau must be of length k=\min(M,N). The matrix Q is related to these components by, Q = Q_k ... Q_2 Q_1 where Q_i = I - \tau_i v_i v_i^T and v_i is the Householder vector

v_i = (0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))

This is the same storage scheme as used by LAPACK. The vector norm is a workspace of length N used for column pivoting.

The algorithm used to perform the decomposition is Householder QR with column pivoting (Golub & Van Loan, “Matrix Computations”, Algorithm 5.4.1).

int gsl_linalg_QRPT_decomp2(const gsl_matrix *A, gsl_matrix *q, gsl_matrix *r, gsl_vector *tau, gsl_permutation *p, int *signum, gsl_vector *norm)

This function factorizes the matrix A into the decomposition A = Q R P^T without modifying A itself and storing the output in the separate matrices q and r.

int gsl_linalg_QRPT_solve(const gsl_matrix *QR, const gsl_vector *tau, const gsl_permutation *p, const gsl_vector *b, gsl_vector *x)

This function solves the square system A x = b using the QRP^T decomposition of A held in (QR, tau, p) which must have been computed previously by gsl_linalg_QRPT_decomp().

int gsl_linalg_QRPT_svx(const gsl_matrix *QR, const gsl_vector *tau, const gsl_permutation *p, gsl_vector *x)

This function solves the square system A x = b in-place using the QRP^T decomposition of A held in (QR, tau, p). On input x should contain the right-hand side b, which is replaced by the solution on output.

int gsl_linalg_QRPT_lssolve(const gsl_matrix *QR, const gsl_vector *tau, const gsl_permutation *p, const gsl_vector *b, gsl_vector *x, gsl_vector *residual)

This function finds the least squares solution to the overdetermined system A x = b where the matrix A has more rows than columns and is assumed to have full rank. The least squares solution minimizes the Euclidean norm of the residual, ||b - A x||. The routine requires as input the QR decomposition of A into (QR, tau, p) given by gsl_linalg_QRPT_decomp(). The solution is returned in x. The residual is computed as a by-product and stored in residual. For rank deficient matrices, gsl_linalg_QRPT_lssolve2() should be used instead.

int gsl_linalg_QRPT_lssolve2(const gsl_matrix *QR, const gsl_vector *tau, const gsl_permutation *p, const gsl_vector *b, const size_t rank, gsl_vector *x, gsl_vector *residual)

This function finds the least squares solution to the overdetermined system A x = b where the matrix A has more rows than columns and has rank given by the input rank. If the user does not know the rank of A, the routine gsl_linalg_QRPT_rank() can be called to estimate it. The least squares solution is the so-called “basic” solution discussed above and may not be the minimum norm solution. The routine requires as input the QR decomposition of A into (QR, tau, p) given by gsl_linalg_QRPT_decomp(). The solution is returned in x. The residual is computed as a by-product and stored in residual.

int gsl_linalg_QRPT_QRsolve(const gsl_matrix *Q, const gsl_matrix *R, const gsl_permutation *p, const gsl_vector *b, gsl_vector *x)

This function solves the square system R P^T x = Q^T b for x. It can be used when the QR decomposition of a matrix is available in unpacked form as (Q, R).

int gsl_linalg_QRPT_update(gsl_matrix *Q, gsl_matrix *R, const gsl_permutation *p, gsl_vector *w, const gsl_vector *v)

This function performs a rank-1 update w v^T of the QRP^T decomposition (Q, R, p). The update is given by Q'R' = Q (R + w v^T P) where the output matrices Q' and R' are also orthogonal and right triangular. Note that w is destroyed by the update. The permutation p is not changed.

int gsl_linalg_QRPT_Rsolve(const gsl_matrix *QR, const gsl_permutation *p, const gsl_vector *b, gsl_vector *x)

This function solves the triangular system R P^T x = b for the N-by-N matrix R contained in QR.

int gsl_linalg_QRPT_Rsvx(const gsl_matrix *QR, const gsl_permutation *p, gsl_vector *x)

This function solves the triangular system R P^T x = b in-place for the N-by-N matrix R contained in QR. On input x should contain the right-hand side b, which is replaced by the solution on output.

size_t gsl_linalg_QRPT_rank(const gsl_matrix *QR, const double tol)

This function estimates the rank of the triangular matrix R contained in QR. The algorithm simply counts the number of diagonal elements of R whose absolute value is greater than the specified tolerance tol. If the input tol is negative, a default value of 20 (M + N) eps(max(|diag(R)|)) is used.

int gsl_linalg_QRPT_rcond(const gsl_matrix *QR, double *rcond, gsl_vector *work)

This function estimates the reciprocal condition number (using the 1-norm) of the R factor, stored in the upper triangle of QR. The reciprocal condition number estimate, defined as 1 / (||R||_1 \cdot ||R^{-1}||_1), is stored in rcond. Additional workspace of size 3 N is required in work.

LQ Decomposition

A general rectangular M-by-N matrix A has a LQ decomposition into the product of a lower trapezoidal M-by-N matrix L and an orthogonal N-by-N square matrix Q:

A = L Q

If M \le N, then L can be written as L = (L_1 \quad 0) where L_1 is M-by-M lower triangular, and

A =
\begin{pmatrix}
  L_1 & 0
\end{pmatrix}
\begin{pmatrix}
  Q_1 \\ Q_2
\end{pmatrix} = L_1 Q_1

where Q_1 consists of the first M rows of Q, and Q_2 contains the remaining N - M rows. The LQ factorization of A is essentially the same as the QR factorization of A^T.

The LQ factorization may be used to find the minimum norm solution of an underdetermined system of equations A x = b, where A is M-by-N and M \le N. The solution is

x = Q^T
\begin{pmatrix}
  L_1^{-1} b \\
  0
\end{pmatrix}

int gsl_linalg_LQ_decomp(gsl_matrix *A, gsl_vector *tau)

This function factorizes the M-by-N matrix A into the LQ decomposition A = L Q. On output the diagonal and lower trapezoidal part of the input matrix contain the matrix L. The vector tau and the elements above the diagonal of the matrix A contain the Householder coefficients and Householder vectors which encode the orthogonal matrix Q. The vector tau must be of length k=\min(M,N). The matrix Q is related to these components by, Q = Q_k ... Q_2 Q_1 where Q_i = I - \tau_i v_i v_i^T and v_i is the Householder vector v_i = (0,...,1,A(i,i+1),A(i,i+2),...,A(i,N)). This is the same storage scheme as used by LAPACK.

int gsl_linalg_LQ_lssolve(const gsl_matrix *LQ, const gsl_vector *tau, const gsl_vector *b, gsl_vector *x, gsl_vector *residual)

This function finds the minimum norm least squares solution to the underdetermined system A x = b, where the M-by-N matrix A has M \le N. The routine requires as input the LQ decomposition of A into (LQ, tau) given by gsl_linalg_LQ_decomp(). The solution is returned in x. The residual, b - Ax, is computed as a by-product and stored in residual.

int gsl_linalg_LQ_unpack(const gsl_matrix *LQ, const gsl_vector *tau, gsl_matrix *Q, gsl_matrix *L)

This function unpacks the encoded LQ decomposition (LQ, tau) into the matrices Q and L, where Q is N-by-N and L is M-by-N.

int gsl_linalg_LQ_QTvec(const gsl_matrix *LQ, const gsl_vector *tau, gsl_vector *v)

This function applies Q^T to the vector v, storing the result Q^T v in v on output.

QL Decomposition

A general rectangular M-by-N matrix A has a QL decomposition into the product of an orthogonal M-by-M square matrix Q (where Q^T Q = I) and an M-by-N left-triangular matrix L.

When M \ge N, the decomposition is given by

A = Q \begin{pmatrix} 0 \\ L_1 \end{pmatrix}

where L_1 is N-by-N lower triangular. When M \le N, the decomposition is given by

A = Q \begin{pmatrix} L_1 & L_2 \end{pmatrix}

where L_1 is a dense M-by-N-M matrix and L_2 is a lower triangular M-by-M matrix.

int gsl_linalg_QL_decomp(gsl_matrix *A, gsl_vector *tau)

This function factorizes the M-by-N matrix A into the QL decomposition A = Q L. The vector tau must be of length N and contains the Householder coefficients on output. The matrix Q is stored in packed form in A on output, using the same storage scheme as LAPACK.

int gsl_linalg_QL_unpack(const gsl_matrix *QL, const gsl_vector *tau, gsl_matrix *Q, gsl_matrix *L)

This function unpacks the encoded QL decomposition (QL, tau) into the matrices Q and L, where Q is M-by-M and L is M-by-N.

Complete Orthogonal Decomposition

The complete orthogonal decomposition of a M-by-N matrix A is a generalization of the QR decomposition with column pivoting, given by

A P = Q
\left(
\begin{matrix}
  R_{11} & 0 \\
  0 & 0
\end{matrix}
\right) Z^T

where P is a N-by-N permutation matrix, Q is M-by-M orthogonal, R_{11} is r-by-r upper triangular, with r = {\rm rank}(A), and Z is N-by-N orthogonal. If A has full rank, then R_{11} = R, Z = I and this reduces to the QR decomposition with column pivoting.

For a rank deficient least squares problem, \min_x{|| b - Ax||^2}, the solution vector x is not unique. However if we further require that ||x||^2 is minimized, then the complete orthogonal decomposition gives us the ability to compute the unique minimum norm solution, which is given by

x = P Z
\left(
\begin{matrix}
  R_{11}^{-1} c_1 \\
  0
\end{matrix}
\right)

and the vector c_1 is the first r elements of Q^T b.

The COD also enables a straightforward solution of regularized least squares problems in Tikhonov standard form, written as

\min_x ||b - A x||^2 + \lambda^2 ||x||^2

where \lambda > 0 is a regularization parameter which represents a tradeoff between minimizing the residual norm ||b - A x|| and the solution norm ||x||. For this system, the solution is given by

x = P Z
\left(
\begin{matrix}
  y_1 \\
  0
\end{matrix}
\right)

where y_1 is a vector of length r which is found by solving

\left(
\begin{matrix}
  R_{11} \\
  \lambda I_r
\end{matrix}
\right) y_1 =
\left(
\begin{matrix}
  c_1 \\
  0
\end{matrix}
\right)

and c_1 is defined above. The equation above can be solved efficiently for different values of \lambda using QR factorizations of the left hand side matrix.

int gsl_linalg_COD_decomp(gsl_matrix *A, gsl_vector *tau_Q, gsl_vector *tau_Z, gsl_permutation *p, size_t *rank, gsl_vector *work)
int gsl_linalg_COD_decomp_e(gsl_matrix *A, gsl_vector *tau_Q, gsl_vector *tau_Z, gsl_permutation *p, double tol, size_t *rank, gsl_vector *work)

These functions factor the M-by-N matrix A into the decomposition A = Q R Z P^T. The rank of A is computed as the number of diagonal elements of R greater than the tolerance tol and output in rank. If tol is not specified, a default value is used (see gsl_linalg_QRPT_rank()). On output, the permutation matrix P is stored in p. The matrix R_{11} is stored in the upper rank-by-rank block of A. The matrices Q and Z are encoded in packed storage in A on output. The vectors tau_Q and tau_Z contain the Householder scalars corresponding to the matrices Q and Z respectively and must be of length k = \min(M,N). The vector work is additional workspace of length N.

int gsl_linalg_COD_lssolve(const gsl_matrix *QRZT, const gsl_vector *tau_Q, const gsl_vector *tau_Z, const gsl_permutation *p, const size_t rank, const gsl_vector *b, gsl_vector *x, gsl_vector *residual)

This function finds the unique minimum norm least squares solution to the overdetermined system A x = b where the matrix A has more rows than columns. The least squares solution minimizes the Euclidean norm of the residual, ||b - A x|| as well as the norm of the solution ||x||. The routine requires as input the QRZT decomposition of A into (QRZT, tau_Q, tau_Z, p, rank) given by gsl_linalg_COD_decomp(). The solution is returned in x. The residual, b - Ax, is computed as a by-product and stored in residual.

int gsl_linalg_COD_lssolve2(const double lambda, const gsl_matrix *QRZT, const gsl_vector *tau_Q, const gsl_vector *tau_Z, const gsl_permutation *p, const size_t rank, const gsl_vector *b, gsl_vector *x, gsl_vector *residual, gsl_matrix *S, gsl_vector *work)

This function finds the solution to the regularized least squares problem in Tikhonov standard form, \min_x ||b - Ax||^2 + \lambda^2 ||x||^2. The routine requires as input the QRZT decomposition of A into (QRZT, tau_Q, tau_Z, p, rank) given by gsl_linalg_COD_decomp(). The parameter \lambda is supplied in lambda. The solution is returned in x. The residual, b - Ax, is stored in residual on output. S is additional workspace of size rank-by-rank. work is additional workspace of length rank.

int gsl_linalg_COD_unpack(const gsl_matrix *QRZT, const gsl_vector *tau_Q, const gsl_vector *tau_Z, const size_t rank, gsl_matrix *Q, gsl_matrix *R, gsl_matrix *Z)

This function unpacks the encoded QRZT decomposition (QRZT, tau_Q, tau_Z, rank) into the matrices Q, R, and Z, where Q is M-by-M, R is M-by-N, and Z is N-by-N.

int gsl_linalg_COD_matZ(const gsl_matrix *QRZT, const gsl_vector *tau_Z, const size_t rank, gsl_matrix *A, gsl_vector *work)

This function multiplies the input matrix A on the right by Z, A' = A Z using the encoded QRZT decomposition (QRZT, tau_Z, rank). A must have N columns but may have any number of rows. Additional workspace of length M is provided in work.

Singular Value Decomposition

A general rectangular M-by-N matrix A has a singular value decomposition (SVD) into the product of an M-by-N orthogonal matrix U, an N-by-N diagonal matrix of singular values S and the transpose of an N-by-N orthogonal square matrix V,

A = U S V^T

The singular values \sigma_i = S_{ii} are all non-negative and are generally chosen to form a non-increasing sequence

\sigma_1 \ge \sigma_2 \ge ... \ge \sigma_N \ge 0

The singular value decomposition of a matrix has many practical uses. The condition number of the matrix is given by the ratio of the largest singular value to the smallest singular value. The presence of a zero singular value indicates that the matrix is singular. The number of non-zero singular values indicates the rank of the matrix. In practice singular value decomposition of a rank-deficient matrix will not produce exact zeroes for singular values, due to finite numerical precision. Small singular values should be edited by choosing a suitable tolerance.

For a rank-deficient matrix, the null space of A is given by the columns of V corresponding to the zero singular values. Similarly, the range of A is given by columns of U corresponding to the non-zero singular values.

Note that the routines here compute the “thin” version of the SVD with U as M-by-N orthogonal matrix. This allows in-place computation and is the most commonly-used form in practice. Mathematically, the “full” SVD is defined with U as an M-by-M orthogonal matrix and S as an M-by-N diagonal matrix (with additional rows of zeros).

int gsl_linalg_SV_decomp(gsl_matrix *A, gsl_matrix *V, gsl_vector *S, gsl_vector *work)

This function factorizes the M-by-N matrix A into the singular value decomposition A = U S V^T for M \ge N. On output the matrix A is replaced by U. The diagonal elements of the singular value matrix S are stored in the vector S. The singular values are non-negative and form a non-increasing sequence from S_1 to S_N. The matrix V contains the elements of V in untransposed form. To form the product U S V^T it is necessary to take the transpose of V. A workspace of length N is required in work.

This routine uses the Golub-Reinsch SVD algorithm.

int gsl_linalg_SV_decomp_mod(gsl_matrix *A, gsl_matrix *X, gsl_matrix *V, gsl_vector *S, gsl_vector *work)

This function computes the SVD using the modified Golub-Reinsch algorithm, which is faster for M \gg N. It requires the vector work of length N and the N-by-N matrix X as additional working space.

int gsl_linalg_SV_decomp_jacobi(gsl_matrix *A, gsl_matrix *V, gsl_vector *S)

This function computes the SVD of the M-by-N matrix A using one-sided Jacobi orthogonalization for M \ge N. The Jacobi method can compute singular values to higher relative accuracy than Golub-Reinsch algorithms (see references for details).

int gsl_linalg_SV_solve(const gsl_matrix *U, const gsl_matrix *V, const gsl_vector *S, const gsl_vector *b, gsl_vector *x)

This function solves the system A x = b using the singular value decomposition (U, S, V) of A which must have been computed previously with gsl_linalg_SV_decomp().

Only non-zero singular values are used in computing the solution. The parts of the solution corresponding to singular values of zero are ignored. Other singular values can be edited out by setting them to zero before calling this function.

In the over-determined case where A has more rows than columns the system is solved in the least squares sense, returning the solution x which minimizes ||A x - b||_2.

int gsl_linalg_SV_solve2(const double tol, const gsl_matrix *U, const gsl_matrix *V, const gsl_vector *S, const gsl_vector *b, gsl_vector *x, gsl_vector *work)

This function solves the system A x = b using the singular value decomposition (U, S, V) of A which must have been computed previously with gsl_linalg_SV_decomp().

Singular values which satisfy, s_i \leq tol \times s_{max} are excluded from the solution. Additional workspace of length N must be provided in work.

In the over-determined case where A has more rows than columns the system is solved in the least squares sense, returning the solution x which minimizes ||b - A x||_2.

int gsl_linalg_SV_lssolve(const double lambda, const gsl_matrix *U, const gsl_matrix *V, const gsl_vector *S, const gsl_vector *b, gsl_vector *x, double *rnorm, gsl_vector *work)

This function solves the regularized least squares problem,

\min_x || b - A x ||^2 + \lambda^2 ||x||^2

using the singular value decomposition (U, S, V) of A which must have been computed previously with gsl_linalg_SV_decomp(). The residual norm ||b - Ax|| is stored in rnorm on output. Additional workspace of size M+N is required in work.

int gsl_linalg_SV_leverage(const gsl_matrix *U, gsl_vector *h)

This function computes the statistical leverage values h_i of a matrix A using its singular value decomposition (U, S, V) previously computed with gsl_linalg_SV_decomp(). h_i are the diagonal values of the matrix A (A^T A)^{-1} A^T and depend only on the matrix U which is the input to this function.

Cholesky Decomposition

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,

A = L L^T

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.

If the matrix A is near singular, it is sometimes possible to reduce the condition number and recover a more accurate solution vector x by scaling as

\left( S A S \right) \left( S^{-1} x \right) = S b

where S is a diagonal matrix whose elements are given by S_{ii} = 1/\sqrt{A_{ii}}. This scaling is also known as Jacobi preconditioning. There are routines below to solve both the scaled and unscaled systems.

int gsl_linalg_cholesky_decomp1(gsl_matrix *A)
int gsl_linalg_complex_cholesky_decomp(gsl_matrix_complex *A)

These functions factorize the symmetric, positive-definite square matrix A into the Cholesky decomposition A = L L^T (or A = L L^{\dagger} for the complex case). On input, the values from the diagonal and lower-triangular part of the matrix A are used (the upper triangular part is ignored). On output the diagonal and lower triangular part of the input matrix A contain the matrix L, while the upper triangular part contains the original matrix. If the matrix is not positive-definite then the decomposition will fail, returning the error code GSL_EDOM.

When testing whether a matrix is positive-definite, disable the error handler first to avoid triggering an error. These functions use Level 3 BLAS to compute the Cholesky factorization (Peise and Bientinesi, 2016).

int gsl_linalg_cholesky_decomp(gsl_matrix *A)

This function is now deprecated and is provided only for backward compatibility.

int gsl_linalg_cholesky_solve(const gsl_matrix *cholesky, const gsl_vector *b, gsl_vector *x)
int gsl_linalg_complex_cholesky_solve(const gsl_matrix_complex *cholesky, const gsl_vector_complex *b, gsl_vector_complex *x)

These functions solve the system A x = b using the Cholesky decomposition of A held in the matrix cholesky which must have been previously computed by gsl_linalg_cholesky_decomp() or gsl_linalg_complex_cholesky_decomp().

int gsl_linalg_cholesky_svx(const gsl_matrix *cholesky, gsl_vector *x)
int gsl_linalg_complex_cholesky_svx(const gsl_matrix_complex *cholesky, gsl_vector_complex *x)

These functions solve the system A x = b in-place using the Cholesky decomposition of A held in the matrix cholesky which must have been previously computed by gsl_linalg_cholesky_decomp() or gsl_linalg_complex_cholesky_decomp(). On input x should contain the right-hand side b, which is replaced by the solution on output.

int gsl_linalg_cholesky_invert(gsl_matrix *cholesky)
int gsl_linalg_complex_cholesky_invert(gsl_matrix_complex *cholesky)

These functions compute the inverse of a matrix from its Cholesky decomposition cholesky, which must have been previously computed by gsl_linalg_cholesky_decomp() or gsl_linalg_complex_cholesky_decomp(). On output, the inverse is stored in-place in cholesky.

int gsl_linalg_cholesky_decomp2(gsl_matrix *A, gsl_vector *S)
int gsl_linalg_complex_cholesky_decomp2(gsl_matrix_complex *A, gsl_vector *S)

This function calculates a diagonal scaling transformation S for the symmetric, positive-definite square matrix A, and then computes the Cholesky decomposition S A S = L L^T. On input, the values from the diagonal and lower-triangular part of the matrix A are used (the upper triangular part is ignored). On output the diagonal and lower triangular part of the input matrix A contain the matrix L. If the matrix is not positive-definite then the decomposition will fail, returning the error code GSL_EDOM. The diagonal scale factors are stored in S on output.

When testing whether a matrix is positive-definite, disable the error handler first to avoid triggering an error.

int gsl_linalg_cholesky_solve2(const gsl_matrix *cholesky, const gsl_vector *S, const gsl_vector *b, gsl_vector *x)
int gsl_linalg_complex_cholesky_solve2(const gsl_matrix_complex *cholesky, const gsl_vector *S, const gsl_vector_complex *b, gsl_vector_complex *x)

This function solves the system (S A S) (S^{-1} x) = S b using the Cholesky decomposition of S A S held in the matrix cholesky which must have been previously computed by gsl_linalg_cholesky_decomp2() or gsl_linalg_complex_cholesky_decomp2().

int gsl_linalg_cholesky_svx2(const gsl_matrix *cholesky, const gsl_vector *S, gsl_vector *x)
int gsl_linalg_complex_cholesky_svx2(const gsl_matrix_complex *cholesky, const gsl_vector *S, gsl_vector_complex *x)

This function solves the system (S A S) (S^{-1} x) = S b in-place using the Cholesky decomposition of S A S held in the matrix cholesky which must have been previously computed by gsl_linalg_cholesky_decomp2() or gsl_linalg_complex_cholesky_decomp2(). On input x should contain the right-hand side b, which is replaced by the solution on output.

int gsl_linalg_cholesky_scale(const gsl_matrix *A, gsl_vector *S)
int gsl_linalg_complex_cholesky_scale(const gsl_matrix_complex *A, gsl_vector *S)

This function calculates a diagonal scaling transformation of the symmetric, positive definite matrix A, such that S A S has a condition number within a factor of N of the matrix of smallest possible condition number over all possible diagonal scalings. On output, S contains the scale factors, given by S_i = 1/\sqrt{A_{ii}}. For any A_{ii} \le 0, the corresponding scale factor S_i is set to 1.

int gsl_linalg_cholesky_scale_apply(gsl_matrix *A, const gsl_vector *S)
int gsl_linalg_complex_cholesky_scale_apply(gsl_matrix_complex *A, const gsl_vector *S)

This function applies the scaling transformation S to the matrix A. On output, A is replaced by S A S.

int gsl_linalg_cholesky_rcond(const gsl_matrix *cholesky, double *rcond, gsl_vector *work)

This function estimates the reciprocal condition number (using the 1-norm) of the symmetric positive definite matrix A, using its Cholesky decomposition provided in cholesky. The reciprocal condition number estimate, defined as 1 / (||A||_1 \cdot ||A^{-1}||_1), is stored in rcond. Additional workspace of size 3 N is required in work.

Pivoted Cholesky Decomposition

A symmetric positive semi-definite square matrix A has an alternate Cholesky decomposition into a product of a lower unit triangular matrix L, a diagonal matrix D and L^T, given by L D L^T. For postive definite matrices, this is equivalent to the Cholesky formulation discussed above, with the standard Cholesky lower triangular factor given by L D^{1 \over 2}. For ill-conditioned matrices, it can help to use a pivoting strategy to prevent the entries of D and L from growing too large, and also ensure D_1 \ge D_2 \ge \cdots \ge D_n > 0, where D_i are the diagonal entries of D. The final decomposition is given by

P A P^T = L D L^T

where P is a permutation matrix.

int gsl_linalg_pcholesky_decomp(gsl_matrix *A, gsl_permutation *p)

This function factors the symmetric, positive-definite square matrix A into the Pivoted Cholesky decomposition P A P^T = L D L^T. On input, the values from the diagonal and lower-triangular part of the matrix A are used to construct the factorization. On output the diagonal of the input matrix A stores the diagonal elements of D, and the lower triangular portion of A contains the matrix L. Since L has ones on its diagonal these do not need to be explicitely stored. The upper triangular portion of A is unmodified. The permutation matrix P is stored in p on output.

int gsl_linalg_pcholesky_solve(const gsl_matrix *LDLT, const gsl_permutation *p, const gsl_vector *b, gsl_vector *x)

This function solves the system A x = b using the Pivoted Cholesky decomposition of A held in the matrix LDLT and permutation p which must have been previously computed by gsl_linalg_pcholesky_decomp().

int gsl_linalg_pcholesky_svx(const gsl_matrix *LDLT, const gsl_permutation *p, gsl_vector *x)

This function solves the system A x = b in-place using the Pivoted Cholesky decomposition of A held in the matrix LDLT and permutation p which must have been previously computed by gsl_linalg_pcholesky_decomp(). On input, x contains the right hand side vector b which is replaced by the solution vector on output.

int gsl_linalg_pcholesky_decomp2(gsl_matrix *A, gsl_permutation *p, gsl_vector *S)

This function computes the pivoted Cholesky factorization of the matrix S A S, where the input matrix A is symmetric and positive definite, and the diagonal scaling matrix S is computed to reduce the condition number of A as much as possible. See Cholesky Decomposition for more information on the matrix S. The Pivoted Cholesky decomposition satisfies P S A S P^T = L D L^T. On input, the values from the diagonal and lower-triangular part of the matrix A are used to construct the factorization. On output the diagonal of the input matrix A stores the diagonal elements of D, and the lower triangular portion of A contains the matrix L. Since L has ones on its diagonal these do not need to be explicitely stored. The upper triangular portion of A is unmodified. The permutation matrix P is stored in p on output. The diagonal scaling transformation is stored in S on output.

int gsl_linalg_pcholesky_solve2(const gsl_matrix *LDLT, const gsl_permutation *p, const gsl_vector *S, const gsl_vector *b, gsl_vector *x)

This function solves the system (S A S) (S^{-1} x) = S b using the Pivoted Cholesky decomposition of S A S held in the matrix LDLT, permutation p, and vector S, which must have been previously computed by gsl_linalg_pcholesky_decomp2().

int gsl_linalg_pcholesky_svx2(const gsl_matrix *LDLT, const gsl_permutation *p, const gsl_vector *S, gsl_vector *x)

This function solves the system (S A S) (S^{-1} x) = S b in-place using the Pivoted Cholesky decomposition of S A S held in the matrix LDLT, permutation p and vector S, which must have been previously computed by gsl_linalg_pcholesky_decomp2(). On input, x contains the right hand side vector b which is replaced by the solution vector on output.

int gsl_linalg_pcholesky_invert(const gsl_matrix *LDLT, const gsl_permutation *p, gsl_matrix *Ainv)

This function computes the inverse of the matrix A, using the Pivoted Cholesky decomposition stored in LDLT and p. On output, the matrix Ainv contains A^{-1}.

int gsl_linalg_pcholesky_rcond(const gsl_matrix *LDLT, const gsl_permutation *p, double *rcond, gsl_vector *work)

This function estimates the reciprocal condition number (using the 1-norm) of the symmetric positive definite matrix A, using its pivoted Cholesky decomposition provided in LDLT. The reciprocal condition number estimate, defined as 1 / (||A||_1 \cdot ||A^{-1}||_1), is stored in rcond. Additional workspace of size 3 N is required in work.

Modified Cholesky Decomposition

The modified Cholesky decomposition is suitable for solving systems A x = b where A is a symmetric indefinite matrix. Such matrices arise in nonlinear optimization algorithms. The standard Cholesky decomposition requires a positive definite matrix and would fail in this case. Instead of resorting to a method like QR or SVD, which do not take into account the symmetry of the matrix, we can instead introduce a small perturbation to the matrix A to make it positive definite, and then use a Cholesky decomposition on the perturbed matrix. The resulting decomposition satisfies

P (A + E) P^T = L D L^T

where P is a permutation matrix, E is a diagonal perturbation matrix, L is unit lower triangular, and D is diagonal. If A is sufficiently positive definite, then the perturbation matrix E will be zero and this method is equivalent to the pivoted Cholesky algorithm. For indefinite matrices, the perturbation matrix E is computed to ensure that A + E is positive definite and well conditioned.

int gsl_linalg_mcholesky_decomp(gsl_matrix *A, gsl_permutation *p, gsl_vector *E)

This function factors the symmetric, indefinite square matrix A into the Modified Cholesky decomposition P (A + E) P^T = L D L^T. On input, the values from the diagonal and lower-triangular part of the matrix A are used to construct the factorization. On output the diagonal of the input matrix A stores the diagonal elements of D, and the lower triangular portion of A contains the matrix L. Since L has ones on its diagonal these do not need to be explicitely stored. The upper triangular portion of A is unmodified. The permutation matrix P is stored in p on output. The diagonal perturbation matrix is stored in E on output. The parameter E may be set to NULL if it is not required.

int gsl_linalg_mcholesky_solve(const gsl_matrix *LDLT, const gsl_permutation *p, const gsl_vector *b, gsl_vector *x)

This function solves the perturbed system (A + E) x = b using the Cholesky decomposition of A + E held in the matrix LDLT and permutation p which must have been previously computed by gsl_linalg_mcholesky_decomp().

int gsl_linalg_mcholesky_svx(const gsl_matrix *LDLT, const gsl_permutation *p, gsl_vector *x)

This function solves the perturbed system (A + E) x = b in-place using the Cholesky decomposition of A + E held in the matrix LDLT and permutation p which must have been previously computed by gsl_linalg_mcholesky_decomp(). On input, x contains the right hand side vector b which is replaced by the solution vector on output.

int gsl_linalg_mcholesky_rcond(const gsl_matrix *LDLT, const gsl_permutation *p, double *rcond, gsl_vector *work)

This function estimates the reciprocal condition number (using the 1-norm) of the perturbed matrix A + E, using its pivoted Cholesky decomposition provided in LDLT. The reciprocal condition number estimate, defined as 1 / (||A + E||_1 \cdot ||(A + E)^{-1}||_1), is stored in rcond. Additional workspace of size 3 N is required in work.

LDLT Decomposition

If A is a symmetric, nonsingular square matrix, then it has a unique factorization of the form

A = L D L^T

where L is a unit lower triangular matrix and D is diagonal. If A is positive definite, then this factorization is equivalent to the Cholesky factorization, where the lower triangular Cholesky factor is L D^{\frac{1}{2}}. Some indefinite matrices for which no Cholesky decomposition exists have an L D L^T decomposition with negative entries in D. The L D L^T algorithm is sometimes referred to as the square root free Cholesky decomposition, as the algorithm does not require the computation of square roots. The algorithm is stable for positive definite matrices, but is not guaranteed to be stable for indefinite matrices.

int gsl_linalg_ldlt_decomp(gsl_matrix *A)

This function factorizes the symmetric, non-singular square matrix A into the decomposition A = L D L^T. On input, the values from the diagonal and lower-triangular part of the matrix A are used. The upper triangle of A is used as temporary workspace. On output the diagonal of A contains the matrix D and the lower triangle of A contains the unit lower triangular matrix L. The matrix 1-norm, ||A||_1 is stored in the upper right corner on output, for later use by gsl_linalg_ldlt_rcond().

If the matrix is detected to be singular, the function returns the error code GSL_EDOM.

int gsl_linalg_ldlt_solve(const gsl_matrix *LDLT, const gsl_vector *b, gsl_vector *x)

This function solves the system A x = b using the L D L^T decomposition of A held in the matrix LDLT which must have been previously computed by gsl_linalg_ldlt_decomp().

int gsl_linalg_ldlt_svx(const gsl_matrix *LDLT, gsl_vector *x)

This function solves the system A x = b in-place using the L D L^T decomposition of A held in the matrix LDLT which must have been previously computed by gsl_linalg_ldlt_decomp(). On input x should contain the right-hand side b, which is replaced by the solution on output.

int gsl_linalg_ldlt_rcond(const gsl_matrix *LDLT, double *rcond, gsl_vector *work)

This function estimates the reciprocal condition number (using the 1-norm) of the symmetric nonsingular matrix A, using its L D L^T decomposition provided in LDLT. The reciprocal condition number estimate, defined as 1 / (||A||_1 \cdot ||A^{-1}||_1), is stored in rcond. Additional workspace of size 3 N is required in work.

Tridiagonal Decomposition of Real Symmetric Matrices

A symmetric matrix A can be factorized by similarity transformations into the form,

A = Q T Q^T

where Q is an orthogonal matrix and T is a symmetric tridiagonal matrix.

int gsl_linalg_symmtd_decomp(gsl_matrix *A, gsl_vector *tau)

This function factorizes the symmetric square matrix A into the symmetric tridiagonal decomposition Q T Q^T. On output the diagonal and subdiagonal part of the input matrix A contain the tridiagonal matrix T. The remaining lower triangular part of the input matrix contains the Householder vectors which, together with the Householder coefficients tau, encode the orthogonal matrix Q. This storage scheme is the same as used by LAPACK. The upper triangular part of A is not referenced.

int gsl_linalg_symmtd_unpack(const gsl_matrix *A, const gsl_vector *tau, gsl_matrix *Q, gsl_vector *diag, gsl_vector *subdiag)

This function unpacks the encoded symmetric tridiagonal decomposition (A, tau) obtained from gsl_linalg_symmtd_decomp() into the orthogonal matrix Q, the vector of diagonal elements diag and the vector of subdiagonal elements subdiag.

int gsl_linalg_symmtd_unpack_T(const gsl_matrix *A, gsl_vector *diag, gsl_vector *subdiag)

This function unpacks the diagonal and subdiagonal of the encoded symmetric tridiagonal decomposition (A, tau) obtained from gsl_linalg_symmtd_decomp() into the vectors diag and subdiag.

Tridiagonal Decomposition of Hermitian Matrices

A hermitian matrix A can be factorized by similarity transformations into the form,

A = U T U^T

where U is a unitary matrix and T is a real symmetric tridiagonal matrix.

int gsl_linalg_hermtd_decomp(gsl_matrix_complex *A, gsl_vector_complex *tau)

This function factorizes the hermitian matrix A into the symmetric tridiagonal decomposition U T U^T. On output the real parts of the diagonal and subdiagonal part of the input matrix A contain the tridiagonal matrix T. The remaining lower triangular part of the input matrix contains the Householder vectors which, together with the Householder coefficients tau, encode the unitary matrix U. This storage scheme is the same as used by LAPACK. The upper triangular part of A and imaginary parts of the diagonal are not referenced.

int gsl_linalg_hermtd_unpack(const gsl_matrix_complex *A, const gsl_vector_complex *tau, gsl_matrix_complex *U, gsl_vector *diag, gsl_vector *subdiag)

This function unpacks the encoded tridiagonal decomposition (A, tau) obtained from gsl_linalg_hermtd_decomp() into the unitary matrix U, the real vector of diagonal elements diag and the real vector of subdiagonal elements subdiag.

int gsl_linalg_hermtd_unpack_T(const gsl_matrix_complex *A, gsl_vector *diag, gsl_vector *subdiag)

This function unpacks the diagonal and subdiagonal of the encoded tridiagonal decomposition (A, tau) obtained from the gsl_linalg_hermtd_decomp() into the real vectors diag and subdiag.

Hessenberg Decomposition of Real Matrices

A general real matrix A can be decomposed by orthogonal similarity transformations into the form

A = U H U^T

where U is orthogonal and H is an upper Hessenberg matrix, meaning that it has zeros below the first subdiagonal. The Hessenberg reduction is the first step in the Schur decomposition for the nonsymmetric eigenvalue problem, but has applications in other areas as well.

int gsl_linalg_hessenberg_decomp(gsl_matrix *A, gsl_vector *tau)

This function computes the Hessenberg decomposition of the matrix A by applying the similarity transformation H = U^T A U. On output, H is stored in the upper portion of A. The information required to construct the matrix U is stored in the lower triangular portion of A. U is a product of N - 2 Householder matrices. The Householder vectors are stored in the lower portion of A (below the subdiagonal) and the Householder coefficients are stored in the vector tau. tau must be of length N.

int gsl_linalg_hessenberg_unpack(gsl_matrix *H, gsl_vector *tau, gsl_matrix *U)

This function constructs the orthogonal matrix U from the information stored in the Hessenberg matrix H along with the vector tau. H and tau are outputs from gsl_linalg_hessenberg_decomp().

int gsl_linalg_hessenberg_unpack_accum(gsl_matrix *H, gsl_vector *tau, gsl_matrix *V)

This function is similar to gsl_linalg_hessenberg_unpack(), except it accumulates the matrix U into V, so that V' = VU. The matrix V must be initialized prior to calling this function. Setting V to the identity matrix provides the same result as gsl_linalg_hessenberg_unpack(). If H is order N, then V must have N columns but may have any number of rows.

int gsl_linalg_hessenberg_set_zero(gsl_matrix *H)

This function sets the lower triangular portion of H, below the subdiagonal, to zero. It is useful for clearing out the Householder vectors after calling gsl_linalg_hessenberg_decomp().

Hessenberg-Triangular Decomposition of Real Matrices

A general real matrix pair (A, B) can be decomposed by orthogonal similarity transformations into the form

A &= U H V^T \\
B &= U R V^T

where U and V are orthogonal, H is an upper Hessenberg matrix, and R is upper triangular. The Hessenberg-Triangular reduction is the first step in the generalized Schur decomposition for the generalized eigenvalue problem.

int gsl_linalg_hesstri_decomp(gsl_matrix *A, gsl_matrix *B, gsl_matrix *U, gsl_matrix *V, gsl_vector *work)

This function computes the Hessenberg-Triangular decomposition of the matrix pair (A, B). On output, H is stored in A, and R is stored in B. If U and V are provided (they may be null), the similarity transformations are stored in them. Additional workspace of length N is needed in work.

Bidiagonalization

A general matrix A can be factorized by similarity transformations into the form,

A = U B V^T

where U and V are orthogonal matrices and B is a N-by-N bidiagonal matrix with non-zero entries only on the diagonal and superdiagonal. The size of U is M-by-N and the size of V is N-by-N.

int gsl_linalg_bidiag_decomp(gsl_matrix *A, gsl_vector *tau_U, gsl_vector *tau_V)

This function factorizes the M-by-N matrix A into bidiagonal form U B V^T. The diagonal and superdiagonal of the matrix B are stored in the diagonal and superdiagonal of A. The orthogonal matrices U and V are stored as compressed Householder vectors in the remaining elements of A. The Householder coefficients are stored in the vectors tau_U and tau_V. The length of tau_U must equal the number of elements in the diagonal of A and the length of tau_V should be one element shorter.

int gsl_linalg_bidiag_unpack(const gsl_matrix *A, const gsl_vector *tau_U, gsl_matrix *U, const gsl_vector *tau_V, gsl_matrix *V, gsl_vector *diag, gsl_vector *superdiag)

This function unpacks the bidiagonal decomposition of A produced by gsl_linalg_bidiag_decomp(), (A, tau_U, tau_V) into the separate orthogonal matrices U, V and the diagonal vector diag and superdiagonal superdiag. Note that U is stored as a compact M-by-N orthogonal matrix satisfying U^T U = I for efficiency.

int gsl_linalg_bidiag_unpack2(gsl_matrix *A, gsl_vector *tau_U, gsl_vector *tau_V, gsl_matrix *V)

This function unpacks the bidiagonal decomposition of A produced by gsl_linalg_bidiag_decomp(), (A, tau_U, tau_V) into the separate orthogonal matrices U, V and the diagonal vector diag and superdiagonal superdiag. The matrix U is stored in-place in A.

int gsl_linalg_bidiag_unpack_B(const gsl_matrix *A, gsl_vector *diag, gsl_vector *superdiag)

This function unpacks the diagonal and superdiagonal of the bidiagonal decomposition of A from gsl_linalg_bidiag_decomp(), into the diagonal vector diag and superdiagonal vector superdiag.

Givens Rotations

A Givens rotation is a rotation in the plane acting on two elements of a given vector. It can be represented in matrix form as

G(i,j,\theta) =
\left(
\begin{matrix}
  1 & \ldots & 0 & \ldots & 0 & \ldots & 0 \\
  \vdots & \ddots & \vdots &  & \vdots & & \vdots \\
  0 & \ldots & \cos{\theta} & \ldots & -\sin{\theta} & \ldots & 0 \\
  \vdots &  & \vdots & \ddots & \vdots & & \vdots \\
  0 & \ldots & \sin{\theta} & \ldots & \cos{\theta} & \ldots & 0 \\
  \vdots &  & \vdots &  & \vdots & \ddots & \vdots \\
  0 & \ldots & 0 & \ldots & 0 & \ldots & 1
\end{matrix}
\right)

where the \cos{\theta} and \sin{\theta} appear at the intersection of the i-th and j-th rows and columns. When acting on a vector x, G(i,j,\theta) x performs a rotation of the (i,j) elements of x. Givens rotations are typically used to introduce zeros in vectors, such as during the QR decomposition of a matrix. In this case, it is typically desired to find c and s such that

\left(
\begin{matrix}
  c & -s \\
  s & c
\end{matrix}
\right)
\left(
\begin{matrix}
  a \\
  b
\end{matrix}
\right) =
\left(
\begin{matrix}
  r \\
  0
\end{matrix}
\right)

with r = \sqrt{a^2 + b^2}.

void gsl_linalg_givens(const double a, const double b, double *c, double *s)

This function computes c = \cos{\theta} and s = \sin{\theta} so that the Givens matrix G(\theta) acting on the vector (a,b) produces (r, 0), with r = \sqrt{a^2 + b^2}.

void gsl_linalg_givens_gv(gsl_vector *v, const size_t i, const size_t j, const double c, const double s)

This function applies the Givens rotation defined by c = \cos{\theta} and s = \sin{\theta} to the i and j elements of v. On output, (v(i),v(j)) \leftarrow G(\theta) (v(i),v(j)).

Householder Transformations

A Householder transformation is a rank-1 modification of the identity matrix which can be used to zero out selected elements of a vector. A Householder matrix H takes the form,

H = I - \tau v v^T

where v is a vector (called the Householder vector) and \tau = 2/(v^T v). The functions described in this section use the rank-1 structure of the Householder matrix to create and apply Householder transformations efficiently.

double gsl_linalg_householder_transform(gsl_vector *w)
gsl_complex gsl_linalg_complex_householder_transform(gsl_vector_complex *w)

This function prepares a Householder transformation H = I - \tau v v^T which can be used to zero all the elements of the input vector w except the first. On output the Householder vector v is stored in w and the scalar \tau is returned. The householder vector v is normalized so that v[0] = 1, however this 1 is not stored in the output vector. Instead, w[0] is set to the first element of the transformed vector, so that if u = H w, w[0] = u[0] on output and the remainder of u is zero.

int gsl_linalg_householder_hm(double tau, const gsl_vector *v, gsl_matrix *A)
int gsl_linalg_complex_householder_hm(gsl_complex tau, const gsl_vector_complex *v, gsl_matrix_complex *A)

This function applies the Householder matrix H defined by the scalar tau and the vector v to the left-hand side of the matrix A. On output the result H A is stored in A.

int gsl_linalg_householder_mh(double tau, const gsl_vector *v, gsl_matrix *A)
int gsl_linalg_complex_householder_mh(gsl_complex tau, const gsl_vector_complex *v, gsl_matrix_complex *A)

This function applies the Householder matrix H defined by the scalar tau and the vector v to the right-hand side of the matrix A. On output the result A H is stored in A.

int gsl_linalg_householder_hv(double tau, const gsl_vector *v, gsl_vector *w)
int gsl_linalg_complex_householder_hv(gsl_complex tau, const gsl_vector_complex *v, gsl_vector_complex *w)

This function applies the Householder transformation H defined by the scalar tau and the vector v to the vector w. On output the result H w is stored in w.

Householder solver for linear systems

int gsl_linalg_HH_solve(gsl_matrix *A, const gsl_vector *b, gsl_vector *x)

This function solves the system A x = b directly using Householder transformations. On output the solution is stored in x and b is not modified. The matrix A is destroyed by the Householder transformations.

int gsl_linalg_HH_svx(gsl_matrix *A, gsl_vector *x)

This function solves the system A x = b in-place using Householder transformations. On input x should contain the right-hand side b, which is replaced by the solution on output. The matrix A is destroyed by the Householder transformations.

Tridiagonal Systems

The functions described in this section efficiently solve symmetric, non-symmetric and cyclic tridiagonal systems with minimal storage. Note that the current implementations of these functions use a variant of Cholesky decomposition, so the tridiagonal matrix must be positive definite. For non-positive definite matrices, the functions return the error code GSL_ESING.

int gsl_linalg_solve_tridiag(const gsl_vector *diag, const gsl_vector *e, const gsl_vector *f, const gsl_vector *b, gsl_vector *x)

This function solves the general N-by-N system A x = b where A is tridiagonal (N \geq 2). The super-diagonal and sub-diagonal vectors e and f must be one element shorter than the diagonal vector diag. The form of A for the 4-by-4 case is shown below,

A =
\left(
\begin{matrix}
  d_0&e_0&  0& 0\\
  f_0&d_1&e_1& 0\\
  0  &f_1&d_2&e_2\\
  0  &0  &f_2&d_3
\end{matrix}
\right)

int gsl_linalg_solve_symm_tridiag(const gsl_vector *diag, const gsl_vector *e, const gsl_vector *b, gsl_vector *x)

This function solves the general N-by-N system A x = b where A is symmetric tridiagonal (N \geq 2). The off-diagonal vector e must be one element shorter than the diagonal vector diag. The form of A for the 4-by-4 case is shown below,

A =
\left(
\begin{matrix}
  d_0&e_0&  0& 0\\
  e_0&d_1&e_1& 0\\
  0  &e_1&d_2&e_2\\
  0  &0  &e_2&d_3
\end{matrix}
\right)

int gsl_linalg_solve_cyc_tridiag(const gsl_vector *diag, const gsl_vector *e, const gsl_vector *f, const gsl_vector *b, gsl_vector *x)

This function solves the general N-by-N system A x = b where A is cyclic tridiagonal (N \geq 3). The cyclic super-diagonal and sub-diagonal vectors e and f must have the same number of elements as the diagonal vector diag. The form of A for the 4-by-4 case is shown below,

A =
\left(
\begin{matrix}
  d_0&e_0& 0 &f_3\\
  f_0&d_1&e_1& 0 \\
  0 &f_1&d_2&e_2\\
  e_3& 0 &f_2&d_3
\end{matrix}
\right)

int gsl_linalg_solve_symm_cyc_tridiag(const gsl_vector *diag, const gsl_vector *e, const gsl_vector *b, gsl_vector *x)

This function solves the general N-by-N system A x = b where A is symmetric cyclic tridiagonal (N \geq 3). The cyclic off-diagonal vector e must have the same number of elements as the diagonal vector diag. The form of A for the 4-by-4 case is shown below,

A =
\left(
\begin{matrix}
  d_0&e_0& 0 &e_3\\
  e_0&d_1&e_1& 0 \\
  0 &e_1&d_2&e_2\\
  e_3& 0 &e_2&d_3
\end{matrix}
\right)

Triangular Systems

int gsl_linalg_tri_invert(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, gsl_matrix *T)
int gsl_linalg_complex_tri_invert(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, gsl_matrix_complex *T)

These functions compute the in-place inverse of the triangular matrix T, stored in the lower triangle when Uplo = CblasLower and upper triangle when Uplo = CblasUpper. The parameter Diag = CblasUnit, CblasNonUnit specifies whether the matrix is unit triangular.

int gsl_linalg_tri_LTL(gsl_matrix *L)
int gsl_linalg_complex_tri_LHL(gsl_matrix_complex *L)

These functions compute the product L^T L (or L^{\dagger} L) in-place and stores it in the lower triangle of L on output.

int gsl_linalg_tri_UL(gsl_matrix *LU)
int gsl_linalg_complex_tri_UL(gsl_matrix_complex *LU)

These functions compute the product U L where U is upper triangular and L is unit lower triangular, stored in LU, as computed by gsl_linalg_LU_decomp() or gsl_linalg_complex_LU_decomp(). The product is computed in-place using Level 3 BLAS.

int gsl_linalg_tri_rcond(CBLAS_UPLO_t Uplo, const gsl_matrix *A, double *rcond, gsl_vector *work)

This function estimates the 1-norm reciprocal condition number of the triangular matrix A, using the lower triangle when Uplo is CblasLower and upper triangle when Uplo is CblasUpper. The reciprocal condition number 1 / \left( \left|\left| A \right|\right|_1 \left|\left| A^{-1} \right|\right|_1 \right) is stored in rcond on output. Additional workspace of size 3N is required in work.

Banded Systems

Band matrices are sparse matrices whose non-zero entries are confined to a diagonal band. From a storage point of view, significant savings can be achieved by storing only the non-zero diagonals of a banded matrix. Algorithms such as LU and Cholesky factorizations preserve the band structure of these matrices. Computationally, working with compact banded matrices is preferable to working on the full dense matrix with many zero entries.

General Banded Format

An example of a general banded matrix is given below.

A = \begin{pmatrix}
      \alpha_1 & \beta_1 & \gamma_1 & 0 & 0 & 0 \\
      \delta_1 & \alpha_2 & \beta_2 & \gamma_2 & 0 & 0 \\
      0 & \delta_2 & \alpha_3 & \beta_3 & \gamma_3 & 0 \\
      0 & 0 & \delta_3 & \alpha_4 & \beta_4 & \gamma_4 \\
      0 & 0 & 0 & \delta_4 & \alpha_5 & \beta_5 \\
      0 & 0 & 0 & 0 & \delta_5 & \alpha_6
    \end{pmatrix}

This matrix has a lower bandwidth of 1 and an upper bandwidth of 2. The lower bandwidth is the number of non-zero subdiagonals, and the upper bandwidth is the number of non-zero superdiagonals. A (p,q) banded matrix has a lower bandwidth p and upper bandwidth q. For example, diagonal matrices are (0,0), tridiagonal matrices are (1,1), and upper triangular matrices are (0,N-1) banded matrices.

The corresponding 6-by-4 packed banded matrix looks like

AB = \begin{pmatrix}
       * & * & \alpha_1 & \delta_1 \\
       * & \beta_1 & \alpha_2 & \delta_2 \\
       \gamma_1 & \beta_2 & \alpha_3 & \delta_3 \\
       \gamma_2 & \beta_3 & \alpha_4 & \delta_4 \\
       \gamma_3 & \beta_4 & \alpha_5 & \delta_5 \\
       \gamma_4 & \beta_5 & \alpha_6 & *
     \end{pmatrix}

where the superdiagonals are stored in columns, followed by the diagonal, followed by the subdiagonals. The entries marked by * are not referenced by the banded routines. With this format, each row of AB corresponds to the non-zero entries of the corresponding column of A. For an N-by-N matrix A, the dimension of AB will be N-by-(p+q+1).

Symmetric Banded Format

Symmetric banded matrices allow for additional storage savings. As an example, consider the following 6 \times 6 symmetric banded matrix with lower bandwidth p = 2:

A = \begin{pmatrix}
      \alpha_1 & \beta_1 & \gamma_1 & 0 & 0 & 0 \\
      \beta_1 & \alpha_2 & \beta_2 & \gamma_2 & 0 & 0 \\
      \gamma_1 & \beta_2 & \alpha_3 & \beta_3 & \gamma_3 & 0 \\
      0 & \gamma_2 & \beta_3 & \alpha_4 & \beta_4 & \gamma_4 \\
      0 & 0 & \gamma_3 & \beta_4 & \alpha_5 & \beta_5 \\
      0 & 0 & 0 & \gamma_4 & \beta_5 & \alpha_6
    \end{pmatrix}

The packed symmetric banded 6 \times 3 matrix will look like:

AB = \begin{pmatrix}
       \alpha_1 & \beta_1 & \gamma_1 \\
       \alpha_2 & \beta_2 & \gamma_2 \\
       \alpha_3 & \beta_3 & \gamma_3 \\
       \alpha_4 & \beta_4 & \gamma_4 \\
       \alpha_5 & \beta_5 & * \\
       \alpha_6 & * & *
     \end{pmatrix}

The entries marked by * are not referenced by the symmetric banded routines. The relationship between the packed format and original matrix is,

AB(i,j) = A(i, i + j) = A(i + j, i)

for i = 0, \dots, N - 1, j = 0, \dots, p. Conversely,

A(i,j) = AB(j, i - j)

for i = 0, \dots, N - 1, j = \textrm{max}(0,i-p), \dots, i.

Warning

Note that this format is the transpose of the symmetric banded format used by LAPACK. In order to develop efficient routines for symmetric banded matrices, it helps to have the nonzero elements in each column in contiguous memory locations. Since C uses row-major order, GSL stores the columns in the rows of the packed banded format, while LAPACK, written in Fortran, uses the transposed format.

Banded LU Decomposition

The routines in this section are designed to factor banded M-by-N matrices with an LU factorization, P A = L U. The matrix A is banded of type (p,q), i.e. a lower bandwidth of p and an upper bandwidth of q. See LU Decomposition for more information on the factorization. For banded (p,q) matrices, the U factor will have an upper bandwidth of p + q, while the L factor will have a lower bandwidth of at most p. Therefore, additional storage is needed to store the p additional bands of U.

As an example, consider the M = N = 7 matrix with lower bandwidth p = 3 and upper bandwidth q = 2,

A = \begin{pmatrix}
      \alpha_1   & \beta_1    & \gamma_1   & 0          & 0          & 0        & 0 \\
      \delta_1   & \alpha_2   & \beta_2    & \gamma_2   & 0          & 0        & 0 \\
      \epsilon_1 & \delta_2   & \alpha_3   & \beta_3    & \gamma_3   & 0        & 0 \\
      \zeta_1    & \epsilon_2 & \delta_3   & \alpha_4   & \beta_4    & \gamma_4 & 0 \\
      0          & \zeta_2    & \epsilon_3 & \delta_4   & \alpha_5   & \beta_5  & \gamma_5 \\
      0          & 0          & \zeta_3    & \epsilon_4 & \delta_5   & \alpha_6 & \beta_6 \\
      0          & 0          & 0          & \zeta_4    & \epsilon_5 & \delta_6 & \alpha_7
    \end{pmatrix}

The corresponding N-by-2p + q + 1 packed banded matrix looks like

AB = \begin{pmatrix}
       * & * & * & * & * & \alpha_1 & \delta_1 & \epsilon_1 & \zeta_1 \\
       * & * & * & * & \beta_1 & \alpha_2 & \delta_2 & \epsilon_2 & \zeta_2 \\
       * & * & * & \gamma_1 & \beta_2 & \alpha_3 & \delta_3 & \epsilon_3 & \zeta_3 \\
       * & * & - & \gamma_2 & \beta_3 & \alpha_4 & \delta_4 & \epsilon_4 & \zeta_4 \\
       * & - & - & \gamma_3 & \beta_4 & \alpha_5 & \delta_5 & \epsilon_5 & * \\
       - & - & - & \gamma_4 & \beta_5 & \alpha_6 & \delta_6 & * & * \\
       \undermat{p}{- & - & -} & \undermat{q}{\gamma_5 & \beta_6} & \alpha_7 & \undermat{p}{* & * & * & }
     \end{pmatrix}

Entries marked with - are used to store the additional p diagonals of the U factor. Entries marked with * are not referenced by the banded routines.

int gsl_linalg_LU_band_decomp(const size_t M, const size_t lb, const size_t ub, gsl_matrix *AB, gsl_vector_uint *piv)

This function computes the LU factorization of the banded matrix AB which is stored in packed band format (see above) and has dimension N-by-2p + q + 1. The number of rows M of the original matrix is provided in M. The lower bandwidth p is provided in lb and the upper bandwidth q is provided in ub. The vector piv has length \textrm{min}(M,N) and stores the pivot indices on output (for 0 \le i < \textrm{min}(M,N), row i of the matrix was interchanged with row piv[i]). On output, AB contains both the L and U factors in packed format.

int gsl_linalg_LU_band_solve(const size_t lb, const size_t ub, const gsl_matrix *LUB, const gsl_vector_uint *piv, const gsl_vector *b, gsl_vector *x)

This function solves the square system Ax = b using the banded LU factorization (LUB, piv) computed by gsl_linalg_LU_band_decomp(). The lower and upper bandwidths are provided in lb and ub respectively. The right hand side vector is provided in b. The solution vector is stored in x on output.

int gsl_linalg_LU_band_svx(const size_t lb, const size_t ub, const gsl_matrix *LUB, const gsl_vector_uint *piv, gsl_vector *x)

This function solves the square system Ax = b in-place, using the banded LU factorization (LUB, piv) computed by gsl_linalg_LU_band_decomp(). The lower and upper bandwidths are provided in lb and ub respectively. On input, the right hand side vector b is provided in x, which is replaced by the solution vector x on output.

int gsl_linalg_LU_band_unpack(const size_t M, const size_t lb, const size_t ub, const gsl_matrix *LUB, const gsl_vector_uint *piv, gsl_matrix *L, gsl_matrix *U)

This function unpacks the banded LU factorization (LUB, piv) previously computed by gsl_linalg_LU_band_decomp() into the matrices L and U. The matrix U has dimension \textrm{min}(M,N)-by-N and stores the upper triangular factor on output. The matrix L has dimension M-by-\textrm{min}(M,N) and stores the matrix P^T L on output.

Banded Cholesky Decomposition

The routines in this section are designed to factor and solve N-by-N linear systems of the form A x = b where A is a banded, symmetric, and positive definite matrix with lower bandwidth p. See Cholesky Decomposition for more information on the factorization. The lower triangular factor of the Cholesky decomposition preserves the same banded structure as the matrix A, enabling an efficient algorithm which overwrites the original matrix with the L factor.

int gsl_linalg_cholesky_band_decomp(gsl_matrix *A)

This function factorizes the symmetric, positive-definite square matrix A into the Cholesky decomposition A = L L^T. The input matrix A is given in symmetric banded format, and has dimensions N-by-(p + 1), where p is the lower bandwidth of the matrix. On output, the entries of A are replaced by the entries of the matrix L in the same format. In addition, the lower right element of A is used to store the matrix 1-norm, used later by gsl_linalg_cholesky_band_rcond() to calculate the reciprocal condition number.

If the matrix is not positive-definite then the decomposition will fail, returning the error code GSL_EDOM. When testing whether a matrix is positive-definite, disable the error handler first to avoid triggering an error.

int gsl_linalg_cholesky_band_solve(const gsl_matrix *LLT, const gsl_vector *b, gsl_vector *x)
int gsl_linalg_cholesky_band_solvem(const gsl_matrix *LLT, const gsl_matrix *B, gsl_matrix *X)

This function solves the symmetric banded system A x = b (or A X = B) using the Cholesky decomposition of A held in the matrix LLT which must have been previously computed by gsl_linalg_cholesky_band_decomp().

int gsl_linalg_cholesky_band_svx(const gsl_matrix *LLT, gsl_vector *x)
int gsl_linalg_cholesky_band_svxm(const gsl_matrix *LLT, gsl_matrix *X)

This function solves the symmetric banded system A x = b (or A X = B) in-place using the Cholesky decomposition of A held in the matrix LLT which must have been previously computed by gsl_linalg_cholesky_band_decomp(). On input x (or X) should contain the right-hand side b (or B), which is replaced by the solution on output.

int gsl_linalg_cholesky_band_invert(const gsl_matrix *LLT, gsl_matrix *Ainv)

This function computes the inverse of a symmetric banded matrix from its Cholesky decomposition LLT, which must have been previously computed by gsl_linalg_cholesky_band_decomp(). On output, the inverse is stored in Ainv, using both the lower and upper portions.

int gsl_linalg_cholesky_band_unpack(const gsl_matrix *LLT, gsl_matrix *L)

This function unpacks the lower triangular Cholesky factor from LLT and stores it in the lower triangular portion of the N-by-N matrix L. The upper triangular portion of L is not referenced.

int gsl_linalg_cholesky_band_scale(const gsl_matrix *A, gsl_vector *S)

This function calculates a diagonal scaling transformation of the symmetric, positive definite banded matrix A, such that S A S has a condition number within a factor of N of the matrix of smallest possible condition number over all possible diagonal scalings. On output, S contains the scale factors, given by S_i = 1/\sqrt{A_{ii}}. For any A_{ii} \le 0, the corresponding scale factor S_i is set to 1.

int gsl_linalg_cholesky_band_scale_apply(gsl_matrix *A, const gsl_vector *S)

This function applies the scaling transformation S to the banded symmetric positive definite matrix A. On output, A is replaced by S A S.

int gsl_linalg_cholesky_band_rcond(const gsl_matrix *LLT, double *rcond, gsl_vector *work)

This function estimates the reciprocal condition number (using the 1-norm) of the symmetric banded positive definite matrix A, using its Cholesky decomposition provided in LLT. The reciprocal condition number estimate, defined as 1 / (||A||_1 \cdot ||A^{-1}||_1), is stored in rcond. Additional workspace of size 3 N is required in work.

Banded LDLT Decomposition

The routines in this section are designed to factor and solve N-by-N linear systems of the form A x = b where A is a banded, symmetric, and non-singular matrix with lower bandwidth p. See LDLT Decomposition for more information on the factorization. The lower triangular factor of the L D L^T decomposition preserves the same banded structure as the matrix A, enabling an efficient algorithm which overwrites the original matrix with the L and D factors.

int gsl_linalg_ldlt_band_decomp(gsl_matrix *A)

This function factorizes the symmetric, non-singular square matrix A into the decomposition A = L D L^T. The input matrix A is given in symmetric banded format, and has dimensions N-by-(p + 1), where p is the lower bandwidth of the matrix. On output, the entries of A are replaced by the entries of the matrices D and L in the same format.

If the matrix is singular then the decomposition will fail, returning the error code GSL_EDOM.

int gsl_linalg_ldlt_band_solve(const gsl_matrix *LDLT, const gsl_vector *b, gsl_vector *x)

This function solves the symmetric banded system A x = b using the L D L^T decomposition of A held in the matrix LDLT which must have been previously computed by gsl_linalg_ldlt_band_decomp().

int gsl_linalg_ldlt_band_svx(const gsl_matrix *LDLT, gsl_vector *x)

This function solves the symmetric banded system A x = b in-place using the L D L^T decomposition of A held in the matrix LDLT which must have been previously computed by gsl_linalg_ldlt_band_decomp(). On input x should contain the right-hand side b, which is replaced by the solution on output.

int gsl_linalg_ldlt_band_unpack(const gsl_matrix *LDLT, gsl_matrix *L, gsl_vector *D)

This function unpacks the unit lower triangular factor L from LDLT and stores it in the lower triangular portion of the N-by-N matrix L. The upper triangular portion of L is not referenced. The diagonal matrix D is stored in the vector D.

int gsl_linalg_ldlt_band_rcond(const gsl_matrix *LDLT, double *rcond, gsl_vector *work)

This function estimates the reciprocal condition number (using the 1-norm) of the symmetric banded nonsingular matrix A, using its L D L^T decomposition provided in LDLT. The reciprocal condition number estimate, defined as 1 / (||A||_1 \cdot ||A^{-1}||_1), is stored in rcond. Additional workspace of size 3 N is required in work.

Balancing

The process of balancing a matrix applies similarity transformations to make the rows and columns have comparable norms. This is useful, for example, to reduce roundoff errors in the solution of eigenvalue problems. Balancing a matrix A consists of replacing A with a similar matrix

A' = D^{-1} A D

where D is a diagonal matrix whose entries are powers of the floating point radix.

int gsl_linalg_balance_matrix(gsl_matrix *A, gsl_vector *D)

This function replaces the matrix A with its balanced counterpart and stores the diagonal elements of the similarity transformation into the vector D.

Examples

The following program solves the linear system A x = b. The system to be solved is,

\left(
\begin{matrix}
  0.18& 0.60& 0.57& 0.96\\
  0.41& 0.24& 0.99& 0.58\\
  0.14& 0.30& 0.97& 0.66\\
  0.51& 0.13& 0.19& 0.85
\end{matrix}
\right)
\left(
\begin{matrix}
  x_0\\
  x_1\\
  x_2\\
  x_3
\end{matrix}
\right)
=
\left(
\begin{matrix}
  1.0\\
  2.0\\
  3.0\\
  4.0
\end{matrix}
\right)

and the solution is found using LU decomposition of the matrix A.

#include <stdio.h>
#include <gsl/gsl_linalg.h>

int
main (void)
{
  double a_data[] = { 0.18, 0.60, 0.57, 0.96,
                      0.41, 0.24, 0.99, 0.58,
                      0.14, 0.30, 0.97, 0.66,
                      0.51, 0.13, 0.19, 0.85 };

  double b_data[] = { 1.0, 2.0, 3.0, 4.0 };

  gsl_matrix_view m
    = gsl_matrix_view_array (a_data, 4, 4);

  gsl_vector_view b
    = gsl_vector_view_array (b_data, 4);

  gsl_vector *x = gsl_vector_alloc (4);

  int s;

  gsl_permutation * p = gsl_permutation_alloc (4);

  gsl_linalg_LU_decomp (&m.matrix, p, &s);

  gsl_linalg_LU_solve (&m.matrix, p, &b.vector, x);

  printf ("x = \n");
  gsl_vector_fprintf (stdout, x, "%g");

  gsl_permutation_free (p);
  gsl_vector_free (x);
  return 0;
}

Here is the output from the program,

x =
-4.05205
-12.6056
1.66091
8.69377

This can be verified by multiplying the solution x by the original matrix A using GNU octave,

octave> A = [ 0.18, 0.60, 0.57, 0.96;
              0.41, 0.24, 0.99, 0.58;
              0.14, 0.30, 0.97, 0.66;
              0.51, 0.13, 0.19, 0.85 ];

octave> x = [ -4.05205; -12.6056; 1.66091; 8.69377];

octave> A * x
ans =
  1.0000
  2.0000
  3.0000
  4.0000

This reproduces the original right-hand side vector, b, in accordance with the equation A x = b.

References and Further Reading

Further information on the algorithms described in this section can be found in the following book,

  • G. H. Golub, C. F. Van Loan, “Matrix Computations” (3rd Ed, 1996), Johns Hopkins University Press, ISBN 0-8018-5414-8.

The LAPACK library is described in the following manual,

  • LAPACK Users’ Guide (Third Edition, 1999), Published by SIAM, ISBN 0-89871-447-8

The LAPACK source code can be found at http://www.netlib.org/lapack, along with an online copy of the users guide.

Further information on recursive Level 3 BLAS algorithms may be found in the following paper,

The recursive Level 3 BLAS QR decomposition is described in the following paper,

  • E. Elmroth and F. G. Gustavson, 2000. Applying recursion to serial and parallel QR factorization leads to better performance. IBM Journal of Research and Development, 44(4), pp.605-624.

The Modified Golub-Reinsch algorithm is described in the following paper,

  • T.F. Chan, “An Improved Algorithm for Computing the Singular Value Decomposition”, ACM Transactions on Mathematical Software, 8 (1982), pp 72–83.

The Jacobi algorithm for singular value decomposition is described in the following papers,

  • J.C. Nash, “A one-sided transformation method for the singular value decomposition and algebraic eigenproblem”, Computer Journal, Volume 18, Number 1 (1975), p 74–76

  • J.C. Nash and S. Shlien “Simple algorithms for the partial singular value decomposition”, Computer Journal, Volume 30 (1987), p 268–275.

  • J. Demmel, K. Veselic, “Jacobi’s Method is more accurate than QR”, Lapack Working Note 15 (LAWN-15), October 1989. Available from netlib, http://www.netlib.org/lapack/ in the lawns or lawnspdf directories.

The algorithm for estimating a matrix condition number is described in the following paper,

  • N. J. Higham, “FORTRAN codes for estimating the one-norm of a real or complex matrix, with applications to condition estimation”, ACM Trans. Math. Soft., vol. 14, no. 4, pp. 381-396, December 1988.