Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Improve svdfact(A, B) docs #15103

Merged
merged 4 commits into from
Feb 17, 2016
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
69 changes: 0 additions & 69 deletions base/docs/helpdb/Base.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1373,20 +1373,6 @@ for all `i` and `j`.
"""
mapslices

"""
svdvals(A)

Returns the singular values of `A`.
"""
svdvals(A)

"""
svdvals(A, B)

Return only the singular values from the generalized singular value decomposition of `A` and `B`.
"""
svdvals(A, B)

"""
issocket(path) -> Bool

Expand Down Expand Up @@ -1595,26 +1581,6 @@ of the problem that is solved on each processor.
"""
peakflops

"""
svd(A, [thin=true]) -> U, S, V

Wrapper around `svdfact` extracting all parts the factorization to a tuple. Direct use of
`svdfact` is therefore generally more efficient. Computes the SVD of `A`, returning `U`,
vector `S`, and `V` such that `A == U*diagm(S)*V'`. If `thin` is `true`, an economy mode
decomposition is returned. The default is to produce a thin decomposition.
"""
svd

"""
svd(A, B) -> U, V, Q, D1, D2, R0

Wrapper around `svdfact` extracting all parts the factorization to a tuple. Direct use of
`svdfact` is therefore generally more efficient. The function returns the generalized SVD of
`A` and `B`, returning `U`, `V`, `Q`, `D1`, `D2`, and `R0` such that `A = U*D1*R0*Q'` and `B =
V*D2*R0*Q'`.
"""
svd(A::AbstractMatrix, B::AbstractMatrix)

"""
ones(type, dims)

Expand Down Expand Up @@ -1774,13 +1740,6 @@ Sum absolute values of elements of an array over the given dimensions.
"""
sumabs(A, dims)

"""
svdvals!(A)

Returns the singular values of `A`, while saving space by overwriting the input.
"""
svdvals!

"""
consume(task, values...)

Expand Down Expand Up @@ -3002,25 +2961,6 @@ are itself and `None` (but `T` itself is not `None`).
"""
isleaftype

"""
svdfact(A, [thin=true]) -> SVD

Compute the Singular Value Decomposition (SVD) of `A` and return an `SVD` object. `U`, `S`,
`V` and `Vt` can be obtained from the factorization `F` with `F[:U]`, `F[:S]`, `F[:V]` and
`F[:Vt]`, such that `A = U*diagm(S)*Vt`. If `thin` is `true`, an economy mode decomposition
is returned. The algorithm produces `Vt` and hence `Vt` is more efficient to extract than
`V`. The default is to produce a thin decomposition.
"""
svdfact(A)

"""
svdfact(A, B) -> GeneralizedSVD

Compute the generalized SVD of `A` and `B`, returning a `GeneralizedSVD` Factorization
object `F`, such that `A = F[:U]*F[:D1]*F[:R0]*F[:Q]'` and `B = F[:V]*F[:D2]*F[:R0]*F[:Q]'`.
"""
svdfact(A, B)

"""
string(xs...)

Expand Down Expand Up @@ -4578,15 +4518,6 @@ interrupt safe. Intended to be called using `do` block syntax as follows:
"""
disable_sigint

"""
svdfact!(A, [thin=true]) -> SVD

`svdfact!` is the same as [`svdfact`](:func:`svdfact`), but saves space by overwriting the
input `A`, instead of creating a copy. If `thin` is `true`, an economy mode decomposition is
returned. The default is to produce a thin decomposition.
"""
svdfact!

"""
hist2d(M, e1, e2) -> (edge1, edge2, counts)

Expand Down
92 changes: 91 additions & 1 deletion base/linalg/svd.jl
Original file line number Diff line number Diff line change
Expand Up @@ -18,13 +18,41 @@ function svdfact!{T<:BlasFloat}(A::StridedMatrix{T}; thin::Bool=true)
end
SVD(u,s,vt)
end
function svdfact{T}(A::StridedVecOrMat{T};thin = true)

"""
svdfact(A, [thin=true]) -> SVD

Compute the singular value decomposition (SVD) of `A` and return an `SVD` object.

`U`, `S`, `V` and `Vt` can be obtained from the factorization `F` with `F[:U]`,
`F[:S]`, `F[:V]` and `F[:Vt]`, such that `A = U*diagm(S)*Vt`.
The algorithm produces `Vt` and hence `Vt` is more efficient to extract than `V`.

If `thin=true` (default), a thin SVD is returned. For a ``M \\times N`` matrix
`A`, `U` is ``M \\times M`` for a full SVD (`thin=false`) and
``M \\times \\min(M, N)`` for a thin SVD.
"""
function svdfact{T}(A::StridedVecOrMat{T}; thin::Bool = true)
S = promote_type(Float32, typeof(one(T)/norm(one(T))))
svdfact!(copy_oftype(A, S), thin = thin)
end
svdfact(x::Number; thin::Bool=true) = SVD(x == 0 ? fill(one(x), 1, 1) : fill(x/abs(x), 1, 1), [abs(x)], fill(one(x), 1, 1))
svdfact(x::Integer; thin::Bool=true) = svdfact(float(x), thin=thin)

"""
svd(A, [thin=true]) -> U, S, V

Computes the SVD of `A`, returning `U`, vector `S`, and `V` such that
`A == U*diagm(S)*V'`.

If `thin=true` (default), a thin SVD is returned. For a ``M \\times N`` matrix
`A`, `U` is ``M \\times M`` for a full SVD (`thin=false`) and
``M \\times \\min(M, N)`` for a thin SVD.

`svd` is a wrapper around [`svdfact`](:func:`svdfact(A)`), extracting all parts
of the `SVD` factorization to a tuple. Direct use of `svdfact` is therefore more
efficient.
"""
function svd(A::Union{Number, AbstractArray}; thin::Bool=true)
F = svdfact(A, thin=thin)
F.U, F.S, F.Vt'
Expand All @@ -44,8 +72,19 @@ function getindex(F::SVD, d::Symbol)
end
end

"""
svdvals!(A)

Returns the singular values of `A`, saving space by overwriting the input.
"""
svdvals!{T<:BlasFloat}(A::StridedMatrix{T}) = any([size(A)...].==0) ? zeros(T, 0) : LAPACK.gesdd!('N', A)[2]
svdvals{T<:BlasFloat}(A::AbstractMatrix{T}) = svdvals!(copy(A))

"""
svdvals(A)

Returns the singular values of `A`.
"""
function svdvals{T}(A::AbstractMatrix{T})
S = promote_type(Float32, typeof(one(T)/norm(one(T))))
svdvals!(copy_oftype(A, S))
Expand Down Expand Up @@ -73,6 +112,16 @@ immutable GeneralizedSVD{T,S} <: Factorization{T}
end
GeneralizedSVD{T}(U::AbstractMatrix{T}, V::AbstractMatrix{T}, Q::AbstractMatrix{T}, a::Vector, b::Vector, k::Int, l::Int, R::AbstractMatrix{T}) = GeneralizedSVD{T,typeof(U)}(U, V, Q, a, b, k, l, R)

"""
svdfact!(A, [thin=true]) -> SVD

`svdfact!` is the same as [`svdfact`](:func:`svdfact`), but saves space by
overwriting the input `A`, instead of creating a copy.

If `thin=true` (default), a thin SVD is returned. For a ``M \\times N`` matrix
`A`, `U` is ``M \\times M`` for a full SVD (`thin=false`) and
``M \\times \\min(M, N)`` for a thin SVD.
"""
function svdfact!{T<:BlasFloat}(A::StridedMatrix{T}, B::StridedMatrix{T})
# xggsvd3 replaced xggsvd in LAPACK 3.6.0
if LAPACK.VERSION[] < v"3.6.0"
Expand All @@ -83,11 +132,45 @@ function svdfact!{T<:BlasFloat}(A::StridedMatrix{T}, B::StridedMatrix{T})
GeneralizedSVD(U, V, Q, a, b, Int(k), Int(l), R)
end
svdfact{T<:BlasFloat}(A::StridedMatrix{T}, B::StridedMatrix{T}) = svdfact!(copy(A),copy(B))

"""
svdfact(A, B) -> GeneralizedSVD

Compute the generalized SVD of `A` and `B`, returning a `GeneralizedSVD` factorization
object `F`, such that `A = F[:U]*F[:D1]*F[:R0]*F[:Q]'` and `B = F[:V]*F[:D2]*F[:R0]*F[:Q]'`.

For an M-by-N matrix `A` and P-by-N matrix `B`,

- `F[:U]` is a M-by-M orthogonal matrix,
- `F[:V]` is a P-by-P orthogonal matrix,
- `F[:Q]` is a N-by-N orthogonal matrix,
- `F[:R0]` is a (K+L)-by-N matrix whose rightmost (K+L)-by-(K+L) block is
nonsingular upper block triangular,
- `F[:D1]` is a M-by-(K+L) diagonal matrix with 1s in the first K entries,
- `F[:D2]` is a P-by-(K+L) matrix whose top right L-by-L block is diagonal,

`K+L` is the effective numerical rank of the matrix `[A; B]`.

The entries of `F[:D1]` and `F[:D2]` are related, as explained in the LAPACK
documentation for the
[generalized SVD](http://www.netlib.org/lapack/lug/node36.html) and the
[xGGSVD3](http://www.netlib.org/lapack/explore-html/d6/db3/dggsvd3_8f.html)
routine which is called underneath (in LAPACK 3.6.0 and newer).
"""
function svdfact{TA,TB}(A::StridedMatrix{TA}, B::StridedMatrix{TB})
S = promote_type(Float32, typeof(one(TA)/norm(one(TA))),TB)
return svdfact!(copy_oftype(A, S), copy_oftype(B, S))
end

"""
svd(A, B) -> U, V, Q, D1, D2, R0

Wrapper around [`svdfact`](:func:`svdfact(A, B)`) extracting all parts of the
factorization to a tuple. Direct use of
`svdfact` is therefore generally more efficient. The function returns the generalized SVD of
`A` and `B`, returning `U`, `V`, `Q`, `D1`, `D2`, and `R0` such that `A = U*D1*R0*Q'` and `B =
V*D2*R0*Q'`.
"""
function svd(A::AbstractMatrix, B::AbstractMatrix)
F = svdfact(A, B)
F[:U], F[:V], F[:Q], F[:D1], F[:D2], F[:R0]
Expand Down Expand Up @@ -141,6 +224,13 @@ function svdvals!{T<:BlasFloat}(A::StridedMatrix{T}, B::StridedMatrix{T})
a[1:k + l] ./ b[1:k + l]
end
svdvals{T<:BlasFloat}(A::StridedMatrix{T},B::StridedMatrix{T}) = svdvals!(copy(A),copy(B))

"""
svdvals(A, B)

Return the generalized singular values from the generalized singular value
decomposition of `A` and `B`.
"""
function svdvals{TA,TB}(A::StridedMatrix{TA}, B::StridedMatrix{TB})
S = promote_type(Float32, typeof(one(TA)/norm(one(TA))), TB)
return svdvals!(copy_oftype(A, S), copy_oftype(B, S))
Expand Down
37 changes: 30 additions & 7 deletions doc/stdlib/linalg.rst
Original file line number Diff line number Diff line change
Expand Up @@ -586,19 +586,29 @@ Linear algebra functions in Julia are largely implemented by calling functions f

.. Docstring generated from Julia source

Compute the Singular Value Decomposition (SVD) of ``A`` and return an ``SVD`` object. ``U``\ , ``S``\ , ``V`` and ``Vt`` can be obtained from the factorization ``F`` with ``F[:U]``\ , ``F[:S]``\ , ``F[:V]`` and ``F[:Vt]``\ , such that ``A = U*diagm(S)*Vt``\ . If ``thin`` is ``true``\ , an economy mode decomposition is returned. The algorithm produces ``Vt`` and hence ``Vt`` is more efficient to extract than ``V``\ . The default is to produce a thin decomposition.
Compute the singular value decomposition (SVD) of ``A`` and return an ``SVD`` object.

``U``\ , ``S``\ , ``V`` and ``Vt`` can be obtained from the factorization ``F`` with ``F[:U]``\ , ``F[:S]``\ , ``F[:V]`` and ``F[:Vt]``\ , such that ``A = U*diagm(S)*Vt``\ . The algorithm produces ``Vt`` and hence ``Vt`` is more efficient to extract than ``V``\ .

If ``thin=true`` (default), a thin SVD is returned. For a :math:`M \times N` matrix ``A``\ , ``U`` is :math:`M \times M` for a full SVD (``thin=false``\ ) and :math:`M \times \min(M, N)` for a thin SVD.

.. function:: svdfact!(A, [thin=true]) -> SVD

.. Docstring generated from Julia source

``svdfact!`` is the same as :func:`svdfact`\ , but saves space by overwriting the input ``A``\ , instead of creating a copy. If ``thin`` is ``true``\ , an economy mode decomposition is returned. The default is to produce a thin decomposition.
``svdfact!`` is the same as :func:`svdfact`\ , but saves space by overwriting the input ``A``\ , instead of creating a copy.

If ``thin=true`` (default), a thin SVD is returned. For a :math:`M \times N` matrix ``A``\ , ``U`` is :math:`M \times M` for a full SVD (``thin=false``\ ) and :math:`M \times \min(M, N)` for a thin SVD.

.. function:: svd(A, [thin=true]) -> U, S, V

.. Docstring generated from Julia source

Wrapper around ``svdfact`` extracting all parts the factorization to a tuple. Direct use of ``svdfact`` is therefore generally more efficient. Computes the SVD of ``A``\ , returning ``U``\ , vector ``S``\ , and ``V`` such that ``A == U*diagm(S)*V'``\ . If ``thin`` is ``true``\ , an economy mode decomposition is returned. The default is to produce a thin decomposition.
Computes the SVD of ``A``\ , returning ``U``\ , vector ``S``\ , and ``V`` such that ``A == U*diagm(S)*V'``\ .

If ``thin=true`` (default), a thin SVD is returned. For a :math:`M \times N` matrix ``A``\ , ``U`` is :math:`M \times M` for a full SVD (``thin=false``\ ) and :math:`M \times \min(M, N)` for a thin SVD.

``svd`` is a wrapper around :func:`svdfact(A)`\ , extracting all parts of the ``SVD`` factorization to a tuple. Direct use of ``svdfact`` is therefore more efficient.

.. function:: svdvals(A)

Expand All @@ -610,25 +620,38 @@ Linear algebra functions in Julia are largely implemented by calling functions f

.. Docstring generated from Julia source

Returns the singular values of ``A``\ , while saving space by overwriting the input.
Returns the singular values of ``A``\ , saving space by overwriting the input.

.. function:: svdfact(A, B) -> GeneralizedSVD

.. Docstring generated from Julia source

Compute the generalized SVD of ``A`` and ``B``\ , returning a ``GeneralizedSVD`` Factorization object ``F``\ , such that ``A = F[:U]*F[:D1]*F[:R0]*F[:Q]'`` and ``B = F[:V]*F[:D2]*F[:R0]*F[:Q]'``\ .
Compute the generalized SVD of ``A`` and ``B``\ , returning a ``GeneralizedSVD`` factorization object ``F``\ , such that ``A = F[:U]*F[:D1]*F[:R0]*F[:Q]'`` and ``B = F[:V]*F[:D2]*F[:R0]*F[:Q]'``\ .

For an M-by-N matrix ``A`` and P-by-N matrix ``B``\ ,

* ``F[:U]`` is a M-by-M orthogonal matrix,
* ``F[:V]`` is a P-by-P orthogonal matrix,
* ``F[:Q]`` is a N-by-N orthogonal matrix,
* ``F[:R0]`` is a (K+L)-by-N matrix whose rightmost (K+L)-by-(K+L) block is nonsingular upper block triangular,
* ``F[:D1]`` is a M-by-(K+L) diagonal matrix with 1s in the first K entries,
* ``F[:D2]`` is a P-by-(K+L) matrix whose top right L-by-L block is diagonal,

``K+L`` is the effective numerical rank of the matrix ``[A; B]``\ .

The entries of ``F[:D1]`` and ``F[:D2]`` are related, as explained in the LAPACK documentation for the `generalized SVD <http://www.netlib.org/lapack/lug/node36.html>`_ and the `xGGSVD3 <http://www.netlib.org/lapack/explore-html/d6/db3/dggsvd3_8f.html>`_ routine which is called underneath (in LAPACK 3.6.0 and newer).

.. function:: svd(A, B) -> U, V, Q, D1, D2, R0

.. Docstring generated from Julia source

Wrapper around ``svdfact`` extracting all parts the factorization to a tuple. Direct use of ``svdfact`` is therefore generally more efficient. The function returns the generalized SVD of ``A`` and ``B``\ , returning ``U``\ , ``V``\ , ``Q``\ , ``D1``\ , ``D2``\ , and ``R0`` such that ``A = U*D1*R0*Q'`` and ``B = V*D2*R0*Q'``\ .
Wrapper around :func:`svdfact(A, B)` extracting all parts of the factorization to a tuple. Direct use of ``svdfact`` is therefore generally more efficient. The function returns the generalized SVD of ``A`` and ``B``\ , returning ``U``\ , ``V``\ , ``Q``\ , ``D1``\ , ``D2``\ , and ``R0`` such that ``A = U*D1*R0*Q'`` and ``B = V*D2*R0*Q'``\ .

.. function:: svdvals(A, B)

.. Docstring generated from Julia source

Return only the singular values from the generalized singular value decomposition of ``A`` and ``B``\ .
Return the generalized singular values from the generalized singular value decomposition of ``A`` and ``B``\ .

.. function:: givens{T}(::T, ::T, ::Integer, ::Integer) -> {Givens, T}

Expand Down