Skip to content

Commit 00fe1eb

Browse files
authored
use _copy_oftype (#256)
* use _copy_oftype * Drop Julia v1.5 support * add tests
1 parent 4e78d85 commit 00fe1eb

File tree

5 files changed

+35
-10
lines changed

5 files changed

+35
-10
lines changed

.github/workflows/ci.yml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -10,7 +10,7 @@ jobs:
1010
fail-fast: false
1111
matrix:
1212
version:
13-
- '1.5'
13+
- '1.6'
1414
- '1'
1515
- '^1.8.0-0'
1616
os:

Project.toml

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
name = "BandedMatrices"
22
uuid = "aae01518-5342-5314-be14-df237901396f"
3-
version = "0.16.13"
3+
version = "0.17"
44

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

1212
[compat]
13-
ArrayLayouts = "0.7, 0.8"
14-
FillArrays = "0.11, 0.12, 0.13"
15-
julia = "1.5"
13+
ArrayLayouts = "0.8"
14+
FillArrays = "0.13"
15+
julia = "1.6"
1616

1717
[extras]
1818
Base64 = "2a0f44e3-6c83-55bd-87e4-b1978d98bd5f"

src/BandedMatrices.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -3,7 +3,7 @@ using Base, FillArrays, ArrayLayouts, LinearAlgebra, SparseArrays, Random
33
using LinearAlgebra.LAPACK
44
import Base: axes, axes1, getproperty, iterate, tail
55
import LinearAlgebra: BlasInt, BlasReal, BlasFloat, BlasComplex, axpy!,
6-
copy_oftype, checksquare, adjoint, transpose, AdjOrTrans, HermOrSym,
6+
checksquare, adjoint, transpose, AdjOrTrans, HermOrSym,
77
_chol!, rot180, dot
88
import LinearAlgebra.BLAS: libblas
99
import LinearAlgebra.LAPACK: liblapack, chkuplo, chktrans
@@ -33,7 +33,7 @@ import ArrayLayouts: MemoryLayout, transposelayout, triangulardata,
3333
triangularlayout, MatLdivVec, hermitianlayout, hermitiandata,
3434
materialize!, BlasMatMulMatAdd, BlasMatMulVecAdd, BlasMatLmulVec, BlasMatLdivVec,
3535
colsupport, rowsupport, symmetricuplo, MatMulMatAdd, MatMulVecAdd,
36-
sublayout, sub_materialize, _fill_lmul!,
36+
sublayout, sub_materialize, _fill_lmul!, _copy_oftype,
3737
reflector!, reflectorApply!, _copyto!, checkdimensions,
3838
_qr!, _qr, _lu!, _lu, _factorize, AbstractTridiagonalLayout, TridiagonalLayout,
3939
BidiagonalLayout, bidiagonaluplo, diagonaldata, supdiagonaldata, subdiagonaldata

src/generic/broadcast.jl

Lines changed: 8 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -913,13 +913,19 @@ axpy!(α, A::AbstractBandedMatrix, dest::AbstractMatrix) = banded_axpy!(α, A, d
913913
for op in (:*, :/)
914914
@eval begin
915915
broadcasted(::BandedStyle, ::typeof($op), a::Zeros, b::AbstractArray) = _broadcasted_zeros(a, b)
916-
broadcasted(::BandedStyle, ::typeof($op), a::Ones{T}, b::AbstractArray{V}) where {T,V} = copy_oftype(b, promote_op(*, T, V))
916+
function broadcasted(::BandedStyle, ::typeof($op), a::AbstractArray{T}, b::Ones{V}) where {T,V}
917+
Base.Broadcast.combine_axes(a, b) # dimension check
918+
_copy_oftype(a, promote_op(*, T, V))
919+
end
917920
end
918921
end
919922

920923
for op in (:*, :\)
921924
@eval begin
922925
broadcasted(::BandedStyle, ::typeof($op), a::AbstractArray, b::Zeros) = _broadcasted_zeros(a, b)
923-
broadcasted(::BandedStyle, ::typeof($op), a::AbstractArray{T}, b::Ones{V}) where {T,V} = copy_oftype(a, promote_op(*, T, V))
926+
function broadcasted(::BandedStyle, ::typeof($op), a::Ones{T}, b::AbstractArray{V}) where {T,V}
927+
Base.Broadcast.combine_axes(a, b) # dimension check
928+
_copy_oftype(b, promote_op(*, T, V))
929+
end
924930
end
925931
end

test/test_broadcasting.jl

Lines changed: 20 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
using BandedMatrices, LinearAlgebra, ArrayLayouts, Test
1+
using BandedMatrices, LinearAlgebra, ArrayLayouts, FillArrays, Test
22
import Base: BroadcastStyle
33
import Base.Broadcast: broadcasted
44
import BandedMatrices: BandedStyle, BandedRows
@@ -541,4 +541,23 @@ import BandedMatrices: BandedStyle, BandedRows
541541
@test bandwidths(B) == bandwidths(A)
542542
@test (A .+ A) .+ A == 3A
543543
end
544+
545+
@testset "ones special case" begin
546+
A = brand(5,4,2,1)
547+
@test Ones(5) .* A == A
548+
@test Ones(5,4) .* A == A
549+
@test Ones(5) .\ A == A
550+
@test Ones(5,4) .\ A == A
551+
@test Ones(5) ./ A == inv.(A)
552+
@test_throws DimensionMismatch Ones(3) .* A
553+
@test_throws DimensionMismatch Ones(3,4) .* A
554+
555+
@test A .* Ones(5) == A
556+
@test A .* Ones(5,4) == A
557+
@test A ./ Ones(5) == A
558+
@test A ./ Ones(5,4) == A
559+
@test A .\ Ones(5) == inv.(A)
560+
@test_throws DimensionMismatch A .* Ones(3)
561+
@test_throws DimensionMismatch A .* Ones(3,4)
562+
end
544563
end

0 commit comments

Comments
 (0)