GenericLinearAlgebra.jl
Documentation for GenericLinearAlgebra.jl
LinearAlgebra.svd!
— Functionsvd!(A[, tol, full])::SVD
A generic singular value decomposition (SVD). The implementation only uses Julia functions so the SVD can be computed for any element type provided that the necessary arithmetic operations are supported by the element type.
tol
: The relative tolerance for determining convergence. The default value iseltype(T)
whereT
is the element type of the input matrix bidiagonal (i.e. after converting the matrix to bidiagonal form).full
: If set totrue
then all the left and right singular vectors are returned. If set tofalse
then only the vectors corresponding to the number of rows and columns of the input matrixA
are returned (the default).
Algorithm
...tomorrow
Example
julia> svd(big.([1 2; 3 4]))
SVD{BigFloat, BigFloat, Matrix{BigFloat}}
U factor:
2×2 Matrix{BigFloat}:
-0.404554 0.914514
-0.914514 -0.404554
singular values:
2-element Vector{BigFloat}:
5.464985704219042650451188493284182533042584640492784181017488774646871847029449
0.3659661906262578204229643842614005434788136943931877734325179702209382149672422
Vt factor:
2×2 Matrix{BigFloat}:
-0.576048 -0.817416
-0.817416 0.576048
LinearAlgebra.svdvals!
— Functionsvdvals!(A [, tol])
Generic computation of singular values.
julia> using LinearAlgebra, GenericLinearAlgebra, Quaternions
julia> n = 20;
julia> H = [big(1)/(i + j - 1) for i in 1:n, j in 1:n]; # The Hilbert matrix
julia> Float64(svdvals(H)[end]/svdvals(Float64.(H))[end] - 1) # The relative error of the LAPACK based solution in 64 bit floating point.
-0.9999999999447275
julia> A = qr([Quaternion(randn(4)...) for i in 1:3, j in 1:3]).Q *
Diagonal([3, 2, 1]) *
qr([Quaternion(randn(4)...) for i in 1:3, j in 1:3]).Q'; # A quaternion matrix with the singular value 1, 2, and 3.
julia> svdvals(A) ≈ [3, 2, 1]
true