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

Fix array return type #70

Merged
merged 2 commits into from
Jun 9, 2022
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
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -11,7 +11,7 @@ SparseArrays = "2f01184e-e22b-5df5-ae63-d93ebab69eaf"
StaticArrays = "90137ffa-7385-5640-81b9-e52037218182"

[compat]
StaticArrays = "0.11, 0.12, 1.0"
StaticArrays = "1"
julia = "1"

[extras]
Expand Down
5 changes: 4 additions & 1 deletion src/IMatrix.jl
Original file line number Diff line number Diff line change
Expand Up @@ -15,6 +15,9 @@ size(A::IMatrix{N}, i::Int) where {N} = N
size(A::IMatrix{N}) where {N} = (N, N)
getindex(A::IMatrix{N,T}, i::Integer, j::Integer) where {N,T} = T(i == j)

Base.:(==)(d1::IMatrix{Na}, d2::IMatrix{Nb}) where {Na,Nb} = Na == Nb
Base.isapprox(d1::IMatrix{Na}, d2::IMatrix{Nb}; kwargs...) where {Na,Nb} = Na == Nb

####### sparse matrix ######
nnz(M::IMatrix{N}) where {N} = N
nonzeros(M::IMatrix{N,T}) where {N,T} = ones(T, N)
Expand All @@ -23,4 +26,4 @@ ishermitian(D::IMatrix) = true
isdense(::IMatrix) = false

similar(::IMatrix{N,Tv}, ::Type{T}) where {N,Tv,T} = IMatrix{N,T}()
copyto!(A::IMatrix{N}, B::IMatrix{N}) where {N} = A
copyto!(A::IMatrix{N}, B::IMatrix{N}) where {N} = A
3 changes: 3 additions & 0 deletions src/PermMatrix.jl
Original file line number Diff line number Diff line change
Expand Up @@ -44,6 +44,9 @@ function PermMatrix(
PermMatrix{Tv,Ti,Vv,Vi}(perm, vals)
end

Base.:(==)(d1::PermMatrix, d2::PermMatrix) = SparseMatrixCSC(d1) == SparseMatrixCSC(d2)
Copy link
Member

Choose a reason for hiding this comment

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

Why can't you compare the value directly? Why covert it to CSC?

Copy link
Member Author

@GiggleLiu GiggleLiu Jun 9, 2022

Choose a reason for hiding this comment

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

If there are two entries being zero, then their order are not important. This logic is too complicated.
If we ignore this corner case, some tests can not pass.

Base.isapprox(d1::PermMatrix, d2::PermMatrix; kwargs...) = isapprox(SparseMatrixCSC(d1), SparseMatrixCSC(d2); kwargs...)
Copy link
Member

Choose a reason for hiding this comment

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

Same as above


################# Array Functions ##################

size(M::PermMatrix) = (length(M.perm), length(M.perm))
Expand Down
64 changes: 37 additions & 27 deletions src/arraymath.jl
Original file line number Diff line number Diff line change
Expand Up @@ -37,38 +37,48 @@ import Base: *, /, ==, +, -, ≈
#+(A::PermMatrix, B::PermMatrix) = PermMatrix(A.dv+B.dv, A.ev+B.ev)
#-(A::PermMatrix, B::PermMatrix) = PermMatrix(A.dv-B.dv, A.ev-B.ev)


const IDP = Union{PermMatrix,IMatrix}

for op in [:+, :-, :(==), :≈]
@eval begin
$op(A::IDP, B::SparseMatrixCSC) = $op(SparseMatrixCSC(A), B)
$op(B::SparseMatrixCSC, A::IDP) = $op(B, SparseMatrixCSC(A))

# intra-IDP
$op(A::PermMatrix, B::IDP) = $op(SparseMatrixCSC(A), SparseMatrixCSC(B))
$op(A::IDP, B::PermMatrix) = $op(SparseMatrixCSC(A), SparseMatrixCSC(B))
$op(A::PermMatrix, B::PermMatrix) = $op(SparseMatrixCSC(A), SparseMatrixCSC(B))
end

# intra-ID
if op in [:+, :-]
@eval begin
$op(d1::Diagonal, d2::IMatrix) = Diagonal($op(d1.diag, diag(d2)))
$op(d1::IMatrix, d2::Diagonal) = Diagonal($op(diag(d1), d2.diag))
end
else
for op in [:+, :-]
for MT in [:IMatrix, :PermMatrix]
@eval begin
$op(d1::IMatrix, d2::Diagonal) = $op(diag(d1), d2.diag)
$op(d1::Diagonal, d2::IMatrix) = $op(d1.diag, diag(d2))
$op(d1::IMatrix{Na}, d2::IMatrix{Nb}) where {Na,Nb} = $op(Na, Nb)
# IMatrix, PermMatrix - SparseMatrixCSC
$op(A::$MT, B::SparseMatrixCSC) = $op(SparseMatrixCSC(A), B)
$op(B::SparseMatrixCSC, A::$MT) = $op(B, SparseMatrixCSC(A))
end
end

@eval begin
# IMatrix, PermMatrix - Diagonal
$op(d1::IMatrix, d2::Diagonal) = Diagonal($op(diag(d1), d2.diag))
$op(d1::Diagonal, d2::IMatrix) = Diagonal($op(d1.diag, diag(d2)))
$op(d1::PermMatrix, d2::Diagonal) = $op(SparseMatrixCSC(d1), d2)
$op(d1::Diagonal, d2::PermMatrix) = $op(d1, SparseMatrixCSC(d2))
# PermMatrix - IMatrix
$op(A::PermMatrix, B::IMatrix) = $op(SparseMatrixCSC(A), SparseMatrixCSC(B))
$op(A::IMatrix, B::PermMatrix) = $op(SparseMatrixCSC(A), SparseMatrixCSC(B))
$op(A::PermMatrix, B::PermMatrix) = $op(SparseMatrixCSC(A), SparseMatrixCSC(B))
end
end

# NOTE: promote 2 at least as an integer
# NOTE: promote to integer
+(d1::IMatrix{Na,Ta}, d2::IMatrix{Nb,Tb}) where {Na,Nb,Ta,Tb} =
d1 == d2 ? Diagonal(fill(promote_type(Ta, Tb, Int)(2), Na)) : throw(DimensionMismatch())
-(d1::IMatrix{Na,Ta}, d2::IMatrix{Nb,Tb}) where {Na,Ta,Nb,Tb} =
d1 == d2 ? spzeros(promote_type(Ta, Tb), Na, Na) : throw(DimensionMismatch())

for MT in [:IMatrix, :PermMatrix]
@eval Base.:(==)(A::$MT, B::SparseMatrixCSC) = SparseMatrixCSC(A) == B
@eval Base.:(==)(A::SparseMatrixCSC, B::$MT) = A == SparseMatrixCSC(B)
end
Base.:(==)(d1::IMatrix, d2::Diagonal) = all(isone, d2.diag)
Base.:(==)(d1::Diagonal, d2::IMatrix) = all(isone, d1.diag)
Base.:(==)(d1::PermMatrix, d2::Diagonal) = SparseMatrixCSC(d1) == SparseMatrixCSC(d2)
Base.:(==)(d1::Diagonal, d2::PermMatrix) = SparseMatrixCSC(d1) == SparseMatrixCSC(d2)
Base.:(==)(A::IMatrix, B::PermMatrix) = SparseMatrixCSC(A) == SparseMatrixCSC(B)
Base.:(==)(A::PermMatrix, B::IMatrix) = SparseMatrixCSC(A) == SparseMatrixCSC(B)

for MT in [:IMatrix, :PermMatrix]
@eval Base.isapprox(A::$MT, B::SparseMatrixCSC; kwargs...) = isapprox(SparseMatrixCSC(A), B)
@eval Base.isapprox(A::SparseMatrixCSC, B::$MT; kwargs...) = isapprox(A, SparseMatrixCSC(B))
@eval Base.isapprox(d1::$MT, d2::Diagonal; kwargs...) = isapprox(diag(d1), d2.diag)
@eval Base.isapprox(d1::Diagonal, d2::$MT; kwargs...) = isapprox(d1.diag, diag(d2))
end
Base.isapprox(A::IMatrix, B::PermMatrix; kwargs...) = isapprox(SparseMatrixCSC(A), SparseMatrixCSC(B); kwargs...)
Base.isapprox(A::PermMatrix, B::IMatrix; kwargs...) = isapprox(SparseMatrixCSC(A), SparseMatrixCSC(B); kwargs...)
5 changes: 3 additions & 2 deletions src/linalg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -95,13 +95,14 @@ end
# to matrix
function *(A::PermMatrix, X::AbstractMatrix)
size(X, 1) == size(A, 2) || throw(DimensionMismatch())
return @views A.vals .* X[A.perm, :] # this may be inefficient for sparse CSC matrix.
return A.vals .* view(X,A.perm,:) # this may be inefficient for sparse CSC matrix.
end

function *(X::AbstractMatrix, A::PermMatrix)
mX, nX = size(X)
nX == size(A, 1) || throw(DimensionMismatch())
return @views (transpose(A.vals) .* X)[:, fast_invperm(A.perm)]
perm = fast_invperm(A.perm)
return transpose(view(A.vals, perm)) .* view(X, :, perm)
end

# NOTE: this is just a temperory fix for v0.7. We should overload mul! in
Expand Down
13 changes: 12 additions & 1 deletion test/linalg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -35,6 +35,9 @@ end
for source_ in Any[p1, sp, ds, dv, pm]
for target in Any[p1, sp, ds, dv, pm]
for source in Any[source_, source_', transpose(source_)]
@test (source == target) == (Matrix(source) == Matrix(target))
@test (source == target) ≈ (Matrix(source) ≈ Matrix(target))
# *
lres = source * target
rres = target * source
flres = Matrix(source) * Matrix(target)
Expand All @@ -49,6 +52,14 @@ end
@test lres |> isdense == false
@test lres |> isdense == false
end

# +, -
lres2 = source + target
rres2 = target - source
flres2 = Matrix(source) + Matrix(target)
frres2 = Matrix(target) - Matrix(source)
@test lres2 ≈ flres2
@test rres2 ≈ frres2
end
end
end
Expand Down Expand Up @@ -133,4 +144,4 @@ end
pop = sparse([2, 3], [4, 5], [1, √2], 5, 5)
zop = Diagonal([0, 1 / 2, 1, -1 / 2, 0])
@test (pop*zop)[2, 4] == -0.5
end
end