Skip to content

Commit

Permalink
New copyto! (#163)
Browse files Browse the repository at this point in the history
* Use _copyto!

* v0.15.1

* _qr

* Update bandedqr.jl
  • Loading branch information
dlfivefifty authored Apr 1, 2020
1 parent 2faf902 commit 31d4e67
Show file tree
Hide file tree
Showing 6 changed files with 18 additions and 21 deletions.
4 changes: 2 additions & 2 deletions Project.toml
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
name = "BandedMatrices"
uuid = "aae01518-5342-5314-be14-df237901396f"
version = "0.15"
version = "0.15.1"

[deps]
ArrayLayouts = "4c555306-a7a7-4459-81d9-ec55ddd5c99a"
Expand All @@ -10,7 +10,7 @@ Random = "9a3f8284-a2c9-5f02-9a11-845980a1fd5c"
SparseArrays = "2f01184e-e22b-5df5-ae63-d93ebab69eaf"

[compat]
ArrayLayouts = "0.2"
ArrayLayouts = "0.2.1"
FillArrays = "0.8"
julia = "1"

Expand Down
3 changes: 2 additions & 1 deletion src/BandedMatrices.jl
Original file line number Diff line number Diff line change
Expand Up @@ -34,7 +34,8 @@ import ArrayLayouts: MemoryLayout, transposelayout, triangulardata,
materialize!, BlasMatMulMatAdd, BlasMatMulVecAdd, BlasMatLmulVec, BlasMatLdivVec,
colsupport, rowsupport, symmetricuplo, MatMulMatAdd, MatMulVecAdd,
sublayout, sub_materialize, _fill_lmul!,
reflector!, reflectorApply!
reflector!, reflectorApply!, _copyto!,
_qr!, _qr, _lu!, _lu, _factorize

import FillArrays: AbstractFill, getindex_value

Expand Down
7 changes: 2 additions & 5 deletions src/banded/BandedLU.jl
Original file line number Diff line number Diff line change
Expand Up @@ -88,13 +88,10 @@ function lu!(A::BandedMatrix{T}, pivot::Union{Val{false}, Val{true}} = Val(true)
return BandedLU{T,typeof(A)}(A, ipiv, zero(BlasInt))
end

lu!(A::AbstractBandedMatrix, pivot::Union{Val{false}, Val{true}} = Val(true); check::Bool = true) =
_lu!(::AbstractBandedLayout, axes, A, pivot::Union{Val{false}, Val{true}} = Val(true); check::Bool = true) =
banded_lufact!(A, pivot; check = check)

function lu(A::Union{AbstractBandedMatrix{T}, AbstractBandedMatrix{Complex{T}}, Adjoint{T,<:AbstractBandedMatrix{T}}, Adjoint{Complex{T},<:AbstractBandedMatrix{Complex{T}}},
Transpose{T,<:AbstractBandedMatrix{T}}, Transpose{Complex{T},<:AbstractBandedMatrix{Complex{T}}}},
pivot::Union{Val{false}, Val{true}} = Val(true);
check::Bool = true) where {T<:Real}
function _lu(::AbstractBandedLayout, axes, A, pivot::Union{Val{false}, Val{true}} = Val(true); check::Bool = true) where {T<:Real}
l,u = bandwidths(A)
lu!(BandedMatrix{float(eltype(A))}(A,(l,l+u)), pivot; check = check)
end
Expand Down
20 changes: 9 additions & 11 deletions src/banded/bandedqr.jl
Original file line number Diff line number Diff line change
@@ -1,24 +1,22 @@
qr(A::AbstractBandedMatrix) = banded_qr(A)
qr(A::BandedSubBandedMatrix) = banded_qr(A)

banded_qr(A) = _banded_qr(axes(A), A)
_qr(::AbstractBandedLayout, ax, A) = _banded_qr(ax, A)
_banded_qr(_, A) = qr!(BandedMatrix{float(eltype(A))}(A, (bandwidth(A,1),bandwidth(A,1)+bandwidth(A,2))))



qr(A::Tridiagonal{T}) where T = qr!(BandedMatrix{float(T)}(A, (1,2)))

qr!(A::BandedMatrix) = banded_qr!(A)
qr!(A::BandedSubBandedMatrix) = banded_qr!(A)
_qr!(::AbstractBandedLayout, _, A) = banded_qr!(A)


function _banded_qr!(R::AbstractMatrix{T}, τ) where T
m,n=size(R)
_banded_qr!(R, τ, min(m - 1 + !(T<:Real), n))
end

function _banded_qr!(R::AbstractMatrix, τ, ncols)
D = bandeddata(R)
l,u = bandwidths(R)
ν = l+u+1
m,n=size(R)

for k = 1:min(m - 1 + !(T<:Real), n)
for k = 1:ncols
x = view(D,u+1:min(ν,m-k+u+1), k)
τk = reflector!(x)
τ[k] = τk
Expand Down Expand Up @@ -79,7 +77,7 @@ function banded_qr_lmul!(adjA::Adjoint, B)
@inbounds begin
for j = 1:nB
cs = colsupport(B,j)
for k = max(1,minimum(cs)-l):min(mA,nA,maximum(cs)+l)
for k = max(1,minimum(cs)-l):min(mA,nA,maximum(cs)+l,length(A.τ))
vBj = B[k,j]
for i = k+1:min(k+l,mB)
vBj += conj(D[i-k+u+1,k])*B[i,j]
Expand Down
2 changes: 1 addition & 1 deletion src/banded/linalg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -71,4 +71,4 @@ end
\(A::Adjoint{<:Any,<:BandedLU}, B::Adjoint{<:Any,<:AbstractVecOrMat}) = A \ copy(B)
\(A::Transpose{<:Any,<:BandedLU}, B::Transpose{<:Any,<:AbstractVecOrMat}) = A \ copy(B)

factorize(A::BandedMatrix) = lu(A)
_factorize(::AbstractBandedLayout, _, A) = lu(A)
3 changes: 2 additions & 1 deletion src/generic/broadcast.jl
Original file line number Diff line number Diff line change
Expand Up @@ -52,7 +52,8 @@ copyto!(dest::AbstractArray, bc::Broadcasted{BandedStyle}) =
# copyto!
##

copyto!(dest::AbstractMatrix, src::AbstractBandedMatrix) = (dest .= src)
_copyto!(::AbstractBandedLayout, ::AbstractBandedLayout, dest::AbstractMatrix, src::AbstractMatrix) =
(dest .= src)


function checkbroadcastband(dest, sizesrc, bndssrc)
Expand Down

2 comments on commit 31d4e67

@dlfivefifty
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@JuliaRegistrator
Copy link

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Registration pull request created: JuliaRegistries/General/11952

After the above pull request is merged, it is recommended that a tag is created on this repository for the registered package version.

This will be done automatically if the Julia TagBot GitHub Action is installed, or can be done manually through the github interface, or via:

git tag -a v0.15.1 -m "<description of version>" 31d4e67df1f7ede3e779b479b1884bd632513af8
git push origin v0.15.1

Please sign in to comment.