Skip to content

Commit

Permalink
Merge pull request #17019 from JuliaLang/ksh/docfactor
Browse files Browse the repository at this point in the history
Added more docs for factorize, cleaned up LDLt
  • Loading branch information
kshyatt authored Jun 29, 2016
2 parents 98d49ab + e976c3a commit 844f284
Show file tree
Hide file tree
Showing 5 changed files with 101 additions and 25 deletions.
9 changes: 0 additions & 9 deletions base/docs/helpdb/Base.jl
Original file line number Diff line number Diff line change
Expand Up @@ -7949,15 +7949,6 @@ Return the index of the first element of `A` for which `predicate` returns `true
"""
findfirst

"""
factorize(A)
Compute a convenient factorization (including LU, Cholesky, Bunch-Kaufman, LowerTriangular,
UpperTriangular) of `A`, based upon the type of the input matrix. The return value can then
be reused for efficient solving of multiple systems. For example: `A=factorize(A); x=A\\b; y=A\\C`.
"""
factorize

"""
promote_rule(type1, type2)
Expand Down
34 changes: 34 additions & 0 deletions base/linalg/dense.jl
Original file line number Diff line number Diff line change
Expand Up @@ -362,6 +362,40 @@ function inv{T}(A::StridedMatrix{T})
return convert(typeof(parent(Ai)), Ai)
end

"""
factorize(A)
Compute a convenient factorization of `A`, based upon the type of the input matrix.
`factorize` checks `A` to see if it is symmetric/triangular/etc. if `A` is passed
as a generic matrix. `factorize` checks every element of `A` to verify/rule out
each property. It will short-circuit as soon as it can rule out symmetry/triangular
structure. The return value can be reused for efficient solving of multiple
systems. For example: `A=factorize(A); x=A\\b; y=A\\C`.
| Properties of `A` | type of factorization |
|:---------------------------|:-----------------------------------------------|
| Positive-definite | Cholesky (see [`cholfact`](:func:`cholfact`)) |
| Dense Symmetric/Hermitian | Bunch-Kaufman (see [`bkfact`](:func:`bkfact`)) |
| Sparse Symmetric/Hermitian | LDLt (see [`ldltfact`](:func:`ldltfact`)) |
| Triangular | Triangular |
| Diagonal | Diagonal |
| Bidiagonal | Bidiagonal |
| Tridiagonal | LU (see [`lufact`](:func:`lufact`)) |
| Symmetric real tridiagonal | LDLt (see [`ldltfact`](:func:`ldltfact`)) |
| General square | LU (see [`lufact`](:func:`lufact`)) |
| General non-square | QR (see [`qrfact`](:func:`qrfact`)) |
If `factorize` is called on a Hermitian positive-definite matrix, for instance, then `factorize`
will return a Cholesky factorization.
Example:
```julia
A = diagm(rand(5)) + diagm(rand(4),1); #A is really bidiagonal
factorize(A) #factorize will check to see that A is already factorized
```
This returns a `5×5 Bidiagonal{Float64}`, which can now be passed to other linear algebra functions
(e.g. eigensolvers) which will use specialized methods for `Bidiagonal` types.
"""
function factorize{T}(A::StridedMatrix{T})
m, n = size(A)
if m == n
Expand Down
4 changes: 2 additions & 2 deletions base/linalg/lq.jl
Original file line number Diff line number Diff line change
Expand Up @@ -21,13 +21,13 @@ LQPackedQ{T}(factors::AbstractMatrix{T}, τ::Vector{T}) = LQPackedQ{T,typeof(fac
lqfact!(A) -> LQ
Compute the LQ factorization of `A`, using the input
matrix as a workspace. See also [`lq`](:func:`lq).
matrix as a workspace. See also [`lq`](:func:`lq`).
"""
lqfact!{T<:BlasFloat}(A::StridedMatrix{T}) = LQ(LAPACK.gelqf!(A)...)
"""
lqfact(A) -> LQ
Compute the LQ factorization of `A`. See also [`lq`](:func:`lq).
Compute the LQ factorization of `A`. See also [`lq`](:func:`lq`).
"""
lqfact{T<:BlasFloat}(A::StridedMatrix{T}) = lqfact!(copy(A))
lqfact(x::Number) = lqfact(fill(x,1,1))
Expand Down
26 changes: 21 additions & 5 deletions base/sparse/cholmod.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1257,9 +1257,13 @@ function cholfact!{Tv}(F::Factor{Tv}, A::Sparse{Tv}; shift::Real=0.0)
end

"""
cholfact!(F::Factor, A::Union{SparseMatrixCSC{<:Real},SparseMatrixCSC{Complex{<:Real}},Symmetric{<:Real,SparseMatrixCSC{<:Real,SuiteSparse_long}},Hermitian{Complex{<:Real},SparseMatrixCSC{Complex{<:Real},SuiteSparse_long}}}; shift = 0.0) -> CHOLMOD.Factor
cholfact!(F::Factor, A; shift = 0.0) -> CHOLMOD.Factor
Compute the Cholesky (``LL'``) factorization of `A`, reusing the symbolic factorization `F`.
`A` must be a `SparseMatrixCSC`, `Symmetric{SparseMatrixCSC}`, or
`Hermitian{SparseMatrixCSC}`. Note that even if `A` doesn't
have the type tag, its structure and values must still be
symmetric/Hermitian.
** Note **
Expand Down Expand Up @@ -1294,9 +1298,13 @@ function cholfact(A::Sparse; shift::Real=0.0,
end

"""
cholfact(::Union{SparseMatrixCSC{<:Real},SparseMatrixCSC{Complex{<:Real}},Symmetric{<:Real,SparseMatrixCSC{<:Real,SuiteSparse_long}},Hermitian{Complex{<:Real},SparseMatrixCSC{Complex{<:Real},SuiteSparse_long}}}; shift = 0.0, perm = Int[]) -> CHOLMOD.Factor
cholfact(A; shift = 0.0, perm = Int[]) -> CHOLMOD.Factor
Compute the Cholesky factorization of a sparse positive definite matrix `A`.
`A` must be a `SparseMatrixCSC`, `Symmetric{SparseMatrixCSC}`, or
`Hermitian{SparseMatrixCSC}`. Note that even if `A` doesn't
have the type tag, its structure and values must still be
symmetric/Hermitian.
A fill-reducing permutation is used.
`F = cholfact(A)` is most frequently used to solve systems of equations with `F\\b`,
but also the methods `diag`, `det`, `logdet` are defined for `F`.
Expand Down Expand Up @@ -1347,9 +1355,13 @@ function ldltfact!{Tv}(F::Factor{Tv}, A::Sparse{Tv}; shift::Real=0.0)
end

"""
ldltfact!(F::Factor, A::Union{SparseMatrixCSC{<:Real},SparseMatrixCSC{Complex{<:Real}},Symmetric{<:Real,SparseMatrixCSC{<:Real,SuiteSparse_long}},Hermitian{Complex{<:Real},SparseMatrixCSC{Complex{<:Real},SuiteSparse_long}}}; shift = 0.0) -> CHOLMOD.Factor
ldltfact!(F::Factor, A; shift = 0.0) -> CHOLMOD.Factor
Compute the ``LDL'`` factorization of `A`, reusing the symbolic factorization `F`.
`A` must be a `SparseMatrixCSC`, `Symmetric{SparseMatrixCSC}`, or
`Hermitian{SparseMatrixCSC}`. Note that even if `A` doesn't
have the type tag, its structure and values must still be
symmetric/Hermitian.
** Note **
Expand Down Expand Up @@ -1384,9 +1396,13 @@ function ldltfact(A::Sparse; shift::Real=0.0,
end

"""
ldltfact(::Union{SparseMatrixCSC{<:Real},SparseMatrixCSC{Complex{<:Real}},Symmetric{<:Real,SparseMatrixCSC{<:Real,SuiteSparse_long}},Hermitian{Complex{<:Real},SparseMatrixCSC{Complex{<:Real},SuiteSparse_long}}}; shift = 0.0, perm=Int[]) -> CHOLMOD.Factor
ldltfact(A; shift = 0.0, perm=Int[]) -> CHOLMOD.Factor
Compute the ``LDL'`` factorization of a sparse symmetric or Hermitian matrix.
Compute the ``LDL'`` factorization of a sparse matrix `A`.
`A` must be a `SparseMatrixCSC`, `Symmetric{SparseMatrixCSC}`, or
`Hermitian{SparseMatrixCSC}`. Note that even if `A` doesn't
have the type tag, its structure and values must still be
symmetric/Hermitian.
A fill-reducing permutation is used.
`F = ldltfact(A)` is most frequently used to solve systems of equations `A*x = b` with `F\\b`.
The returned factorization object `F` also supports the methods `diag`,
Expand Down
53 changes: 44 additions & 9 deletions doc/stdlib/linalg.rst
Original file line number Diff line number Diff line change
Expand Up @@ -51,7 +51,42 @@ Linear algebra functions in Julia are largely implemented by calling functions f

.. Docstring generated from Julia source
Compute a convenient factorization (including LU, Cholesky, Bunch-Kaufman, LowerTriangular, UpperTriangular) of ``A``\ , based upon the type of the input matrix. The return value can then be reused for efficient solving of multiple systems. For example: ``A=factorize(A); x=A\b; y=A\C``\ .
Compute a convenient factorization of ``A``\ , based upon the type of the input matrix. ``factorize`` checks ``A`` to see if it is symmetric/triangular/etc. if ``A`` is passed as a generic matrix. ``factorize`` checks every element of ``A`` to verify/rule out each property. It will short-circuit as soon as it can rule out symmetry/triangular structure. The return value can be reused for efficient solving of multiple systems. For example: ``A=factorize(A); x=A\b; y=A\C``\ .

+----------------------------+--------------------------------------+
| Properties of ``A`` | type of factorization |
+============================+======================================+
| Positive-definite | Cholesky (see :func:`cholfact`\ ) |
+----------------------------+--------------------------------------+
| Dense Symmetric/Hermitian | Bunch-Kaufman (see :func:`bkfact`\ ) |
+----------------------------+--------------------------------------+
| Sparse Symmetric/Hermitian | LDLt (see :func:`ldltfact`\ ) |
+----------------------------+--------------------------------------+
| Triangular | Triangular |
+----------------------------+--------------------------------------+
| Diagonal | Diagonal |
+----------------------------+--------------------------------------+
| Bidiagonal | Bidiagonal |
+----------------------------+--------------------------------------+
| Tridiagonal | LU (see :func:`lufact`\ ) |
+----------------------------+--------------------------------------+
| Symmetric real tridiagonal | LDLt (see :func:`ldltfact`\ ) |
+----------------------------+--------------------------------------+
| General square | LU (see :func:`lufact`\ ) |
+----------------------------+--------------------------------------+
| General non-square | QR (see :func:`qrfact`\ ) |
+----------------------------+--------------------------------------+

If ``factorize`` is called on a Hermitian positive-definite matrix, for instance, then ``factorize`` will return a Cholesky factorization.

Example:

.. code-block:: julia
A = diagm(rand(5)) + diagm(rand(4),1); #A is really bidiagonal
factorize(A) #factorize will check to see that A is already factorized
This returns a ``5×5 Bidiagonal{Float64}``\ , which can now be passed to other linear algebra functions (e.g. eigensolvers) which will use specialized methods for ``Bidiagonal`` types.

.. function:: full(F)

Expand Down Expand Up @@ -283,11 +318,11 @@ Linear algebra functions in Julia are largely implemented by calling functions f
Compute the pivoted Cholesky factorization of a dense symmetric positive semi-definite matrix ``A`` and return a ``CholeskyPivoted`` factorization. The ``uplo`` argument may be ``:L`` for using the lower part or ``:U`` for the upper part of ``A``\ . The default is to use ``:U``\ . The triangular Cholesky factor can be obtained from the factorization ``F`` with: ``F[:L]`` and ``F[:U]``\ . The following functions are available for ``PivotedCholesky`` objects: ``size``\ , ``\``\ , ``inv``\ , ``det``\ , and ``rank``\ . The argument ``tol`` determines the tolerance for determining the rank. For negative values, the tolerance is the machine precision.

.. function:: cholfact(::Union{SparseMatrixCSC{<:Real},SparseMatrixCSC{Complex{<:Real}},Symmetric{<:Real,SparseMatrixCSC{<:Real,SuiteSparse_long}},Hermitian{Complex{<:Real},SparseMatrixCSC{Complex{<:Real},SuiteSparse_long}}}; shift = 0.0, perm = Int[]) -> CHOLMOD.Factor
.. function:: cholfact(A; shift = 0.0, perm = Int[]) -> CHOLMOD.Factor

.. Docstring generated from Julia source
Compute the Cholesky factorization of a sparse positive definite matrix ``A``\ . A fill-reducing permutation is used. ``F = cholfact(A)`` is most frequently used to solve systems of equations with ``F\b``\ , but also the methods ``diag``\ , ``det``\ , ``logdet`` are defined for ``F``\ . You can also extract individual factors from ``F``\ , using ``F[:L]``\ . However, since pivoting is on by default, the factorization is internally represented as ``A == P'*L*L'*P`` with a permutation matrix ``P``\ ; using just ``L`` without accounting for ``P`` will give incorrect answers. To include the effects of permutation, it's typically preferable to extact "combined" factors like ``PtL = F[:PtL]`` (the equivalent of ``P'*L``\ ) and ``LtP = F[:UP]`` (the equivalent of ``L'*P``\ ).
Compute the Cholesky factorization of a sparse positive definite matrix ``A``\ . ``A`` must be a ``SparseMatrixCSC``\ , ``Symmetric{SparseMatrixCSC}``\ , or ``Hermitian{SparseMatrixCSC}``\ . A fill-reducing permutation is used. ``F = cholfact(A)`` is most frequently used to solve systems of equations with ``F\b``\ , but also the methods ``diag``\ , ``det``\ , ``logdet`` are defined for ``F``\ . You can also extract individual factors from ``F``\ , using ``F[:L]``\ . However, since pivoting is on by default, the factorization is internally represented as ``A == P'*L*L'*P`` with a permutation matrix ``P``\ ; using just ``L`` without accounting for ``P`` will give incorrect answers. To include the effects of permutation, it's typically preferable to extact "combined" factors like ``PtL = F[:PtL]`` (the equivalent of ``P'*L``\ ) and ``LtP = F[:UP]`` (the equivalent of ``L'*P``\ ).

Setting optional ``shift`` keyword argument computes the factorization of ``A+shift*I`` instead of ``A``\ . If the ``perm`` argument is nonempty, it should be a permutation of ``1:size(A,1)`` giving the ordering to use (instead of CHOLMOD's default AMD ordering).

Expand All @@ -297,11 +332,11 @@ Linear algebra functions in Julia are largely implemented by calling functions f

Many other functions from CHOLMOD are wrapped but not exported from the ``Base.SparseArrays.CHOLMOD`` module.

.. function:: cholfact!(F::Factor, A::Union{SparseMatrixCSC{<:Real},SparseMatrixCSC{Complex{<:Real}},Symmetric{<:Real,SparseMatrixCSC{<:Real,SuiteSparse_long}},Hermitian{Complex{<:Real},SparseMatrixCSC{Complex{<:Real},SuiteSparse_long}}}; shift = 0.0) -> CHOLMOD.Factor
.. function:: cholfact!(F::Factor, A; shift = 0.0) -> CHOLMOD.Factor

.. Docstring generated from Julia source
Compute the Cholesky (:math:`LL'`\ ) factorization of ``A``\ , reusing the symbolic factorization ``F``\ .
Compute the Cholesky (:math:`LL'`\ ) factorization of ``A``\ , reusing the symbolic factorization ``F``\ . ``A`` must be a ``SparseMatrixCSC``\ , ``Symmetric{SparseMatrixCSC}``\ , or ``Hermitian{SparseMatrixCSC}``\ .

** Note **

Expand Down Expand Up @@ -353,11 +388,11 @@ Linear algebra functions in Julia are largely implemented by calling functions f
Compute an ``LDLt`` factorization of a real symmetric tridiagonal matrix such that ``A = L*Diagonal(d)*L'`` where ``L`` is a unit lower triangular matrix and ``d`` is a vector. The main use of an ``LDLt`` factorization ``F = ldltfact(A)`` is to solve the linear system of equations ``Ax = b`` with ``F\b``\ .

.. function:: ldltfact(::Union{SparseMatrixCSC{<:Real},SparseMatrixCSC{Complex{<:Real}},Symmetric{<:Real,SparseMatrixCSC{<:Real,SuiteSparse_long}},Hermitian{Complex{<:Real},SparseMatrixCSC{Complex{<:Real},SuiteSparse_long}}}; shift = 0.0, perm=Int[]) -> CHOLMOD.Factor
.. function:: ldltfact(A; shift = 0.0, perm=Int[]) -> CHOLMOD.Factor

.. Docstring generated from Julia source
Compute the :math:`LDL'` factorization of a sparse symmetric or Hermitian matrix. A fill-reducing permutation is used. ``F = ldltfact(A)`` is most frequently used to solve systems of equations ``A*x = b`` with ``F\b``\ . The returned factorization object ``F`` also supports the methods ``diag``\ , ``det``\ , and ``logdet``\ . You can extract individual factors from ``F`` using ``F[:L]``\ . However, since pivoting is on by default, the factorization is internally represented as ``A == P'*L*D*L'*P`` with a permutation matrix ``P``\ ; using just ``L`` without accounting for ``P`` will give incorrect answers. To include the effects of permutation, it's typically preferable to extact "combined" factors like ``PtL = F[:PtL]`` (the equivalent of ``P'*L``\ ) and ``LtP = F[:UP]`` (the equivalent of ``L'*P``\ ). The complete list of supported factors is ``:L, :PtL, :D, :UP, :U, :LD, :DU, :PtLD, :DUP``\ .
Compute the :math:`LDL'` factorization of a sparse matrix ``A``\ . ``A`` must be a ``SparseMatrixCSC``\ , ``Symmetric{SparseMatrixCSC}``\ , or ``Hermitian{SparseMatrixCSC}``\ . Note that even if ``A`` doesn't have the type tag, its structure must still be symmetric/Hermitian. A fill-reducing permutation is used. ``F = ldltfact(A)`` is most frequently used to solve systems of equations ``A*x = b`` with ``F\b``\ . The returned factorization object ``F`` also supports the methods ``diag``\ , ``det``\ , and ``logdet``\ . You can extract individual factors from ``F`` using ``F[:L]``\ . However, since pivoting is on by default, the factorization is internally represented as ``A == P'*L*D*L'*P`` with a permutation matrix ``P``\ ; using just ``L`` without accounting for ``P`` will give incorrect answers. To include the effects of permutation, it's typically preferable to extact "combined" factors like ``PtL = F[:PtL]`` (the equivalent of ``P'*L``\ ) and ``LtP = F[:UP]`` (the equivalent of ``L'*P``\ ). The complete list of supported factors is ``:L, :PtL, :D, :UP, :U, :LD, :DU, :PtLD, :DUP``\ .
Setting optional ``shift`` keyword argument computes the factorization of ``A+shift*I`` instead of ``A``\ . If the ``perm`` argument is nonempty, it should be a permutation of ``1:size(A,1)`` giving the ordering to use (instead of CHOLMOD's default AMD ordering).

Expand All @@ -367,11 +402,11 @@ Linear algebra functions in Julia are largely implemented by calling functions f

Many other functions from CHOLMOD are wrapped but not exported from the ``Base.SparseArrays.CHOLMOD`` module.

.. function:: ldltfact!(F::Factor, A::Union{SparseMatrixCSC{<:Real},SparseMatrixCSC{Complex{<:Real}},Symmetric{<:Real,SparseMatrixCSC{<:Real,SuiteSparse_long}},Hermitian{Complex{<:Real},SparseMatrixCSC{Complex{<:Real},SuiteSparse_long}}}; shift = 0.0) -> CHOLMOD.Factor
.. function:: ldltfact!(F::Factor, A; shift = 0.0) -> CHOLMOD.Factor

.. Docstring generated from Julia source
Compute the :math:`LDL'` factorization of ``A``\ , reusing the symbolic factorization ``F``\ .
Compute the :math:`LDL'` factorization of ``A``\ , reusing the symbolic factorization ``F``\ . ``A`` must be a ``SparseMatrixCSC``\ , ``Symmetric{SparseMatrixCSC}``\ , or ``Hermitian{SparseMatrixCSC}``\ .

** Note **

Expand Down

0 comments on commit 844f284

Please sign in to comment.