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

Rename rank one up- and downdate functions for Cholesky to rank(up/down)date #14424

Merged
merged 1 commit into from
Dec 18, 2015
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
3 changes: 3 additions & 0 deletions NEWS.md
Original file line number Diff line number Diff line change
Expand Up @@ -70,6 +70,9 @@ Library improvements
appropriate. The `sparsevec` function returns a one-dimensional sparse
vector instead of a one-column sparse matrix. ([#13440])

* Rank one update and downdate functions, `lowrankupdate`, `lowrankupdate!`, `lowrankdowndate`,
and `lowrankdowndate!`, for dense Cholesky factorizations ([#14243],[#14424])

* New `foreach` function for calling a function on every element of a collection when
the results are not needed.

Expand Down
16 changes: 8 additions & 8 deletions base/linalg/cholesky.jl
Original file line number Diff line number Diff line change
Expand Up @@ -267,11 +267,11 @@ chkfullrank(C::CholeskyPivoted) = C.rank < size(C.factors, 1) && throw(RankDefic
rank(C::CholeskyPivoted) = C.rank

"""
update!(C::Cholesky, v::StridedVector) -> CC::Cholesky
lowrankupdate!(C::Cholesky, v::StridedVector) -> CC::Cholesky

Update a Cholesky factorization `C` with the vector `v`. If `A = C[:U]'C[:U]` then `CC = cholfact(C[:U]'C[:U] + v*v')` but the computation of `CC` only uses `O(n^2)` operations. The input factorization `C` is updated in place such that on exit `C == CC`. The vector `v` is destroyed during the computation.
"""
function update!(C::Cholesky, v::StridedVector)
function lowrankupdate!(C::Cholesky, v::StridedVector)
A = C.factors
n = length(v)
if size(C, 1) != n
Expand Down Expand Up @@ -310,11 +310,11 @@ function update!(C::Cholesky, v::StridedVector)
end

"""
downdate!(C::Cholesky, v::StridedVector) -> CC::Cholesky
lowrankdowndate!(C::Cholesky, v::StridedVector) -> CC::Cholesky

Downdate a Cholesky factorization `C` with the vector `v`. If `A = C[:U]'C[:U]` then `CC = cholfact(C[:U]'C[:U] - v*v')` but the computation of `CC` only uses `O(n^2)` operations. The input factorization `C` is updated in place such that on exit `C == CC`. The vector `v` is destroyed during the computation.
"""
function downdate!(C::Cholesky, v::StridedVector)
function lowrankdowndate!(C::Cholesky, v::StridedVector)
A = C.factors
n = length(v)
if size(C, 1) != n
Expand Down Expand Up @@ -360,15 +360,15 @@ function downdate!(C::Cholesky, v::StridedVector)
end

"""
update(C::Cholesky, v::StridedVector) -> CC::Cholesky
lowrankupdate(C::Cholesky, v::StridedVector) -> CC::Cholesky

Update a Cholesky factorization `C` with the vector `v`. If `A = C[:U]'C[:U]` then `CC = cholfact(C[:U]'C[:U] + v*v')` but the computation of `CC` only uses `O(n^2)` operations.
"""
update(C::Cholesky, v::StridedVector) = update!(copy(C), copy(v))
lowrankupdate(C::Cholesky, v::StridedVector) = lowrankupdate!(copy(C), copy(v))

"""
downdate(C::Cholesky, v::StridedVector) -> CC::Cholesky
lowrankdowndate(C::Cholesky, v::StridedVector) -> CC::Cholesky

Downdate a Cholesky factorization `C` with the vector `v`. If `A = C[:U]'C[:U]` then `CC = cholfact(C[:U]'C[:U] - v*v')` but the computation of `CC` only uses `O(n^2)` operations.
"""
downdate(C::Cholesky, v::StridedVector) = downdate!(copy(C), copy(v))
lowrankdowndate(C::Cholesky, v::StridedVector) = lowrankdowndate!(copy(C), copy(v))
8 changes: 4 additions & 4 deletions doc/stdlib/linalg.rst
Original file line number Diff line number Diff line change
Expand Up @@ -161,25 +161,25 @@ Linear algebra functions in Julia are largely implemented by calling functions f

.. currentmodule:: Base.LinAlg

.. function:: update(C::Cholesky, v::StridedVector) -> CC::Cholesky
.. function:: lowrankupdate(C::Cholesky, v::StridedVector) -> CC::Cholesky

.. Docstring generated from Julia source

Update a Cholesky factorization ``C`` with the vector ``v``\ . If ``A = C[:U]'C[:U]`` then ``CC = cholfact(C[:U]'C[:U] + v*v')`` but the computation of ``CC`` only uses ``O(n^2)`` operations.

.. function:: downdate(C::Cholesky, v::StridedVector) -> CC::Cholesky
.. function:: lowrankdowndate(C::Cholesky, v::StridedVector) -> CC::Cholesky

.. Docstring generated from Julia source

Downdate a Cholesky factorization ``C`` with the vector ``v``\ . If ``A = C[:U]'C[:U]`` then ``CC = cholfact(C[:U]'C[:U] - v*v')`` but the computation of ``CC`` only uses ``O(n^2)`` operations.

.. function:: update!(C::Cholesky, v::StridedVector) -> CC::Cholesky
.. function:: lowrankupdate!(C::Cholesky, v::StridedVector) -> CC::Cholesky

.. Docstring generated from Julia source

Update a Cholesky factorization ``C`` with the vector ``v``\ . If ``A = C[:U]'C[:U]`` then ``CC = cholfact(C[:U]'C[:U] + v*v')`` but the computation of ``CC`` only uses ``O(n^2)`` operations. The input factorization ``C`` is updated in place such that on exit ``C == CC``\ . The vector ``v`` is destroyed during the computation.

.. function:: downdate!(C::Cholesky, v::StridedVector) -> CC::Cholesky
.. function:: lowrankdowndate!(C::Cholesky, v::StridedVector) -> CC::Cholesky

.. Docstring generated from Julia source

Expand Down
4 changes: 2 additions & 2 deletions test/linalg/cholesky.jl
Original file line number Diff line number Diff line change
Expand Up @@ -152,7 +152,7 @@ let A = complex(randn(10,5), randn(10, 5)), v = complex(randn(5), randn(5))
for uplo in (:U, :L)
AcA = A'A
F = cholfact(AcA, uplo)
@test LinAlg.update(F, v)[uplo] cholfact(AcA + v*v')[uplo]
@test LinAlg.downdate(F, v)[uplo] cholfact(AcA - v*v')[uplo]
@test LinAlg.lowrankupdate(F, v)[uplo] cholfact(AcA + v*v')[uplo]
@test LinAlg.lowrankdowndate(F, v)[uplo] cholfact(AcA - v*v')[uplo]
end
end