Go to the first, previous, next, last section, table of contents.

Singular Value Decomposition

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

The singular values @math{\sigma_i = S_{ii}} are all non-negative and are generally chosen to form a non-increasing sequence @math{\sigma_1 >= \sigma_2 >= ... >= \sigma_N >= 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.

Function: int gsl_linalg_SV_decomp (gsl_matrix * A, gsl_matrix * V, gsl_vector * S, gsl_vector * work)
This function factorizes the @math{M}-by-@math{N} matrix A into the singular value decomposition @math{A = U S V^T}. On output the matrix A is replaced by @math{U}. The diagonal elements of the singular value matrix @math{S} are stored in the vector S. The singular values are non-negative and form a non-increasing sequence from @math{S_1} to @math{S_N}. The matrix V contains the elements of @math{V} in untransposed form. To form the product @math{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.

Function: 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 @math{M>>N}. It requires the vector work and the @math{N}-by-@math{N} matrix X as additional working space.

Function: int gsl_linalg_SV_decomp_jacobi (gsl_matrix * A, gsl_matrix * V, gsl_vector * S)
This function computes the SVD using one-sided Jacobi orthogonalization (see references for details). The Jacobi method can compute singular values to higher relative accuracy than Golub-Reinsch algorithms.

Function: int gsl_linalg_SV_solve (gsl_matrix * U, gsl_matrix * V, gsl_vector * S, const gsl_vector * b, gsl_vector * x)
This function solves the system @math{A x = b} using the singular value decomposition (U, S, V) of @math{A} given by 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 @math{||A x - b||_2}.


Go to the first, previous, next, last section, table of contents.