Skip to content

Commit 08e0504

Browse files
committed
Use Braille to display adjoint and transposed sparse matrices
1 parent 98e678f commit 08e0504

File tree

2 files changed

+152
-44
lines changed

2 files changed

+152
-44
lines changed

stdlib/SparseArrays/src/sparsematrix.jl

Lines changed: 60 additions & 22 deletions
Original file line numberDiff line numberDiff line change
@@ -110,6 +110,8 @@ julia> nnz(A)
110110
"""
111111
nnz(S::AbstractSparseMatrixCSC) = Int(getcolptr(S)[size(S, 2) + 1] - 1)
112112
nnz(S::ReshapedArray{<:Any,1,<:AbstractSparseMatrixCSC}) = nnz(parent(S))
113+
nnz(S::Adjoint{<:Any,<:AbstractSparseMatrixCSC}) = nnz(parent(S))
114+
nnz(S::Transpose{<:Any,<:AbstractSparseMatrixCSC}) = nnz(parent(S))
113115
nnz(S::UpperTriangular{<:Any,<:AbstractSparseMatrixCSC}) = nnz1(S)
114116
nnz(S::LowerTriangular{<:Any,<:AbstractSparseMatrixCSC}) = nnz1(S)
115117
nnz(S::SparseMatrixCSCView) = nnz1(S)
@@ -200,6 +202,7 @@ nzrange(S::SparseMatrixCSCView, col::Integer) = nzrange(S.parent, S.indices[2][c
200202
nzrange(S::UpperTriangular{<:Any,<:SparseMatrixCSCUnion}, i::Integer) = nzrangeup(S.data, i)
201203
nzrange(S::LowerTriangular{<:Any,<:SparseMatrixCSCUnion}, i::Integer) = nzrangelo(S.data, i)
202204

205+
const AbstractSparseMatrixCSCInclAdjointAndTranspose = Union{AbstractSparseMatrixCSC,Adjoint{<:Any,<:AbstractSparseMatrixCSC},Transpose{<:Any,<:AbstractSparseMatrixCSC}}
203206
function Base.isstored(A::AbstractSparseMatrixCSC, i::Integer, j::Integer)
204207
@boundscheck checkbounds(A, i, j)
205208
rows = rowvals(A)
@@ -209,10 +212,19 @@ function Base.isstored(A::AbstractSparseMatrixCSC, i::Integer, j::Integer)
209212
return false
210213
end
211214

212-
Base.replace_in_print_matrix(A::AbstractSparseMatrix, i::Integer, j::Integer, s::AbstractString) =
215+
function Base.isstored(A::Union{Adjoint{<:Any,<:AbstractSparseMatrixCSC},Transpose{<:Any,<:AbstractSparseMatrixCSC}}, i::Integer, j::Integer)
216+
@boundscheck checkbounds(A, i, j)
217+
cols = rowvals(parent(A))
218+
for istored in nzrange(parent(A), i)
219+
j == cols[istored] && return true
220+
end
221+
return false
222+
end
223+
224+
Base.replace_in_print_matrix(A::AbstractSparseMatrixCSCInclAdjointAndTranspose, i::Integer, j::Integer, s::AbstractString) =
213225
Base.isstored(A, i, j) ? s : Base.replace_with_centered_mark(s)
214226

215-
function Base.show(io::IO, ::MIME"text/plain", S::AbstractSparseMatrixCSC)
227+
function Base.show(io::IO, ::MIME"text/plain", S::AbstractSparseMatrixCSCInclAdjointAndTranspose)
216228
xnnz = nnz(S)
217229
m, n = size(S)
218230
print(io, m, "×", n, " ", typeof(S), " with ", xnnz, " stored ",
@@ -226,7 +238,7 @@ end
226238
Base.show(io::IO, S::AbstractSparseMatrixCSC) = Base.show(convert(IOContext, io), S::AbstractSparseMatrixCSC)
227239

228240
const brailleBlocks = UInt16['', '', '', '', '', '', '', '']
229-
function _show_with_braille_patterns(io::IOContext, S::AbstractSparseMatrixCSC)
241+
function _show_with_braille_patterns(io::IOContext, S::AbstractSparseMatrixCSCInclAdjointAndTranspose)
230242
m, n = size(S)
231243
(m == 0 || n == 0) && return show(io, MIME("text/plain"), S)
232244

@@ -260,28 +272,43 @@ function _show_with_braille_patterns(io::IOContext, S::AbstractSparseMatrixCSC)
260272
brailleGrid = fill(UInt16(10240), (scaleWidth - 1) ÷ 2 + 2, (scaleHeight - 1) ÷ 4 + 1)
261273
brailleGrid[end, :] .= '\n'
262274

263-
rvals = rowvals(S)
275+
rvals = rowvals(parent(S))
264276
rowscale = max(1, scaleHeight - 1) / max(1, m - 1)
265277
colscale = max(1, scaleWidth - 1) / max(1, n - 1)
266-
@inbounds for j = 1:n
267-
# Scale the column index `j` to the best matching column index
268-
# of a matrix of size `scaleHeight × scaleWidth`
269-
sj = round(Int, (j - 1) * colscale + 1)
270-
for x in nzrange(S, j)
271-
# Scale the row index `i` to the best matching row index
278+
if isa(S, AbstractSparseMatrixCSC)
279+
@inbounds for j = 1:n
280+
# Scale the column index `j` to the best matching column index
272281
# of a matrix of size `scaleHeight × scaleWidth`
273-
si = round(Int, (rvals[x] - 1) * rowscale + 1)
274-
275-
# Given the index pair `(si, sj)` of the scaled matrix,
276-
# calculate the corresponding triple `(k, l, p)` such that the
277-
# element at `(si, sj)` can be found at position `(k, l)` in the
278-
# braille grid `brailleGrid` and corresponds to the 1-dot braille
279-
# character `brailleBlocks[p]`
280-
k = (sj - 1) ÷ 2 + 1
281-
l = (si - 1) ÷ 4 + 1
282-
p = ((sj - 1) % 2) * 4 + ((si - 1) % 4 + 1)
283-
284-
brailleGrid[k, l] |= brailleBlocks[p]
282+
sj = round(Int, (j - 1) * colscale + 1)
283+
for x in nzrange(S, j)
284+
# Scale the row index `i` to the best matching row index
285+
# of a matrix of size `scaleHeight × scaleWidth`
286+
si = round(Int, (rvals[x] - 1) * rowscale + 1)
287+
288+
# Given the index pair `(si, sj)` of the scaled matrix,
289+
# calculate the corresponding triple `(k, l, p)` such that the
290+
# element at `(si, sj)` can be found at position `(k, l)` in the
291+
# braille grid `brailleGrid` and corresponds to the 1-dot braille
292+
# character `brailleBlocks[p]`
293+
k = (sj - 1) ÷ 2 + 1
294+
l = (si - 1) ÷ 4 + 1
295+
p = ((sj - 1) % 2) * 4 + ((si - 1) % 4 + 1)
296+
297+
brailleGrid[k, l] |= brailleBlocks[p]
298+
end
299+
end
300+
else
301+
# If `S` is a adjoint or transpose of a sparse matrix we invert the
302+
# roles of the indices `i` and `j`
303+
@inbounds for i = 1:m
304+
si = round(Int, (i - 1) * rowscale + 1)
305+
for x in nzrange(parent(S), i)
306+
sj = round(Int, (rvals[x] - 1) * colscale + 1)
307+
k = (sj - 1) ÷ 2 + 1
308+
l = (si - 1) ÷ 4 + 1
309+
p = ((sj - 1) % 2) * 4 + ((si - 1) % 4 + 1)
310+
brailleGrid[k, l] |= brailleBlocks[p]
311+
end
285312
end
286313
end
287314
foreach(c -> print(io, Char(c)), @view brailleGrid[1:end-1])
@@ -298,6 +325,17 @@ function Base.show(io::IOContext, S::AbstractSparseMatrixCSC)
298325
_show_with_braille_patterns(io, S)
299326
end
300327

328+
function Base.show(io::IOContext, S::Union{Adjoint{<:Any,<:AbstractSparseMatrixCSC},Transpose{<:Any,<:AbstractSparseMatrixCSC}})
329+
if max(size(S)...) < 16 && !(get(io, :compact, false)::Bool)
330+
ioc = IOContext(io, :compact => true)
331+
println(ioc)
332+
Base.print_matrix(ioc, S)
333+
return
334+
end
335+
println(io)
336+
_show_with_braille_patterns(io, S)
337+
end
338+
301339
## Reshape
302340

303341
function sparse_compute_reshaped_colptr_and_rowval(colptrS::Vector{Ti}, rowvalS::Vector{Ti},

stdlib/SparseArrays/test/sparse.jl

Lines changed: 92 additions & 22 deletions
Original file line numberDiff line numberDiff line change
@@ -2281,22 +2281,70 @@ end
22812281
for c in unstored_indices
22822282
@test Base.isstored(A, c[1], c[2]) == false
22832283
end
2284+
2285+
# `isstored` for adjoint and tranposed matrices:
2286+
for trans in (adjoint, transpose)
2287+
B = trans(A)
2288+
stored_indices = [CartesianIndex(j, i) for (j, i) in zip(J, I)]
2289+
unstored_indices = [c for c in CartesianIndices((n, m)) if !(c in stored_indices)]
2290+
for c in stored_indices
2291+
@test Base.isstored(B, c[1], c[2]) == true
2292+
end
2293+
for c in unstored_indices
2294+
@test Base.isstored(B, c[1], c[2]) == false
2295+
end
2296+
end
22842297
end
22852298

22862299
@testset "show" begin
22872300
io = IOBuffer()
2288-
show(io, MIME"text/plain"(), spzeros(Float64, Int64, 0, 0))
2289-
@test String(take!(io)) == "0×0 SparseArrays.SparseMatrixCSC{Float64,Int64} with 0 stored entries"
2290-
show(io, MIME"text/plain"(), sparse(Int64[1], Int64[1], [1.0]))
2291-
@test String(take!(io)) == "1×1 SparseArrays.SparseMatrixCSC{Float64,Int64} with 1 stored entry:\n 1.0"
2292-
show(io, MIME"text/plain"(), spzeros(Float32, Int64, 2, 2))
2293-
@test String(take!(io)) == "2×2 SparseArrays.SparseMatrixCSC{Float32,Int64} with 0 stored entries:\n ⋅ ⋅ \n ⋅ ⋅ "
2301+
2302+
A = spzeros(Float64, Int64, 0, 0)
2303+
for (transform, showstring) in zip(
2304+
(identity, adjoint, transpose), (
2305+
"0×0 SparseArrays.SparseMatrixCSC{Float64,Int64} with 0 stored entries",
2306+
"0×0 LinearAlgebra.Adjoint{Float64,SparseArrays.SparseMatrixCSC{Float64,Int64}} with 0 stored entries",
2307+
"0×0 LinearAlgebra.Transpose{Float64,SparseArrays.SparseMatrixCSC{Float64,Int64}} with 0 stored entries"
2308+
))
2309+
show(io, MIME"text/plain"(), transform(A))
2310+
@test String(take!(io)) == showstring
2311+
end
2312+
2313+
A = sparse(Int64[1], Int64[1], [1.0])
2314+
for (transform, showstring) in zip(
2315+
(identity, adjoint, transpose), (
2316+
"1×1 SparseArrays.SparseMatrixCSC{Float64,Int64} with 1 stored entry:\n 1.0",
2317+
"1×1 LinearAlgebra.Adjoint{Float64,SparseArrays.SparseMatrixCSC{Float64,Int64}} with 1 stored entry:\n 1.0",
2318+
"1×1 LinearAlgebra.Transpose{Float64,SparseArrays.SparseMatrixCSC{Float64,Int64}} with 1 stored entry:\n 1.0",
2319+
))
2320+
show(io, MIME"text/plain"(), transform(A))
2321+
@test String(take!(io)) == showstring
2322+
end
2323+
2324+
A = spzeros(Float32, Int64, 2, 2)
2325+
for (transform, showstring) in zip(
2326+
(identity, adjoint, transpose), (
2327+
"2×2 SparseArrays.SparseMatrixCSC{Float32,Int64} with 0 stored entries:\n ⋅ ⋅ \n ⋅ ⋅ ",
2328+
"2×2 LinearAlgebra.Adjoint{Float32,SparseArrays.SparseMatrixCSC{Float32,Int64}} with 0 stored entries:\n ⋅ ⋅ \n ⋅ ⋅ ",
2329+
"2×2 LinearAlgebra.Transpose{Float32,SparseArrays.SparseMatrixCSC{Float32,Int64}} with 0 stored entries:\n ⋅ ⋅ \n ⋅ ⋅ ",
2330+
))
2331+
show(io, MIME"text/plain"(), transform(A))
2332+
@test String(take!(io)) == showstring
2333+
end
22942334

22952335
A = sparse(Int64[1, 1], Int64[1, 2], [1.0, 2.0])
2296-
show(io, MIME"text/plain"(), A)
2297-
@test String(take!(io)) == "1×2 SparseArrays.SparseMatrixCSC{Float64,Int64} with 2 stored entries:\n 1.0 2.0"
2298-
_show_with_braille_patterns(convert(IOContext, io), A)
2299-
@test String(take!(io)) == ""
2336+
for (transform, showstring, braille) in zip(
2337+
(identity, adjoint, transpose), (
2338+
"1×2 SparseArrays.SparseMatrixCSC{Float64,Int64} with 2 stored entries:\n 1.0 2.0",
2339+
"2×1 LinearAlgebra.Adjoint{Float64,SparseArrays.SparseMatrixCSC{Float64,Int64}} with 2 stored entries:\n 1.0\n 2.0",
2340+
"2×1 LinearAlgebra.Transpose{Float64,SparseArrays.SparseMatrixCSC{Float64,Int64}} with 2 stored entries:\n 1.0\n 2.0",
2341+
),
2342+
("", "", ""))
2343+
show(io, MIME"text/plain"(), transform(A))
2344+
@test String(take!(io)) == showstring
2345+
_show_with_braille_patterns(convert(IOContext, io), transform(A))
2346+
@test String(take!(io)) == braille
2347+
end
23002348

23012349
# every 1-dot braille pattern
23022350
for (i, b) in enumerate(split("⠁⠂⠄⡀⠈⠐⠠⢀", ""))
@@ -2308,25 +2356,47 @@ end
23082356

23092357
# empty braille pattern Char(10240)
23102358
A = spzeros(Int64, Int64, 4, 2)
2311-
_show_with_braille_patterns(convert(IOContext, io), A)
2312-
@test String(take!(io)) == "" * Char(10240)
2359+
for (transform, braille) in zip(
2360+
(identity, adjoint, transpose),
2361+
("" * Char(10240), "" * Char(10240)^2, "" * Char(10240)^2))
2362+
_show_with_braille_patterns(convert(IOContext, io), transform(A))
2363+
@test String(take!(io)) == braille
2364+
end
23132365

23142366
A = sparse(Int64[1, 2, 4, 2, 3], Int64[1, 1, 1, 2, 2], Int64[1, 1, 1, 1, 1], 4, 2)
2315-
show(io, MIME"text/plain"(), A)
2316-
@test String(take!(io)) == "4×2 SparseArrays.SparseMatrixCSC{Int64,Int64} with 5 stored entries:\n 1 ⋅\n 1 1\n ⋅ 1\n 1 ⋅"
2317-
_show_with_braille_patterns(convert(IOContext, io), A)
2318-
@test String(take!(io)) == ""
2367+
for (transform, showstring, braille) in zip(
2368+
(identity, adjoint, transpose), (
2369+
"4×2 SparseArrays.SparseMatrixCSC{Int64,Int64} with 5 stored entries:\n 1 ⋅\n 1 1\n ⋅ 1\n 1 ⋅",
2370+
"2×4 LinearAlgebra.Adjoint{Int64,SparseArrays.SparseMatrixCSC{Int64,Int64}} with 5 stored entries:\n 1 1 ⋅ 1\n ⋅ 1 1 ⋅",
2371+
"2×4 LinearAlgebra.Transpose{Int64,SparseArrays.SparseMatrixCSC{Int64,Int64}} with 5 stored entries:\n 1 1 ⋅ 1\n ⋅ 1 1 ⋅",
2372+
),
2373+
("", "⠙⠊", "⠙⠊"))
2374+
show(io, MIME"text/plain"(), transform(A))
2375+
@test String(take!(io)) == showstring
2376+
_show_with_braille_patterns(convert(IOContext, io), transform(A))
2377+
@test String(take!(io)) == braille
2378+
end
23192379

23202380
A = sparse(Int64[1, 3, 2, 4], Int64[1, 1, 2, 2], Int64[1, 1, 1, 1], 7, 3)
2321-
show(io, MIME"text/plain"(), A)
2322-
@test String(take!(io)) == "7×3 SparseArrays.SparseMatrixCSC{Int64,Int64} with 4 stored entries:\n 1 ⋅ ⋅\n ⋅ 1 ⋅\n 1 ⋅ ⋅\n ⋅ 1 ⋅\n ⋅ ⋅ ⋅\n ⋅ ⋅ ⋅\n ⋅ ⋅ ⋅"
2323-
_show_with_braille_patterns(convert(IOContext, io), A)
2324-
@test String(take!(io)) == "" * Char(10240) * "\n" * Char(10240)^2
2381+
for (transform, showstring, braille) in zip(
2382+
(identity, adjoint, transpose), (
2383+
"7×3 SparseArrays.SparseMatrixCSC{Int64,Int64} with 4 stored entries:\n 1 ⋅ ⋅\n ⋅ 1 ⋅\n 1 ⋅ ⋅\n ⋅ 1 ⋅\n ⋅ ⋅ ⋅\n ⋅ ⋅ ⋅\n ⋅ ⋅ ⋅",
2384+
"3×7 LinearAlgebra.Adjoint{Int64,SparseArrays.SparseMatrixCSC{Int64,Int64}} with 4 stored entries:\n 1 ⋅ 1 ⋅ ⋅ ⋅ ⋅\n ⋅ 1 ⋅ 1 ⋅ ⋅ ⋅\n ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅",
2385+
"3×7 LinearAlgebra.Transpose{Int64,SparseArrays.SparseMatrixCSC{Int64,Int64}} with 4 stored entries:\n 1 ⋅ 1 ⋅ ⋅ ⋅ ⋅\n ⋅ 1 ⋅ 1 ⋅ ⋅ ⋅\n ⋅ ⋅ ⋅ ⋅ ⋅ ⋅ ⋅",
2386+
),
2387+
("" * Char(10240) * "\n" * Char(10240)^2, "⠑⠑" * Char(10240)^2, "⠑⠑" * Char(10240)^2))
2388+
show(io, MIME"text/plain"(), transform(A))
2389+
@test String(take!(io)) == showstring
2390+
_show_with_braille_patterns(convert(IOContext, io), transform(A))
2391+
@test String(take!(io)) == braille
2392+
end
23252393

23262394
A = sparse(Int64[1:10;], Int64[1:10;], fill(Float64(1), 10))
2327-
_show_with_braille_patterns(convert(IOContext, io), A)
23282395
brailleString = "⠑⢄" * Char(10240)^3 * "\n" * Char(10240)^2 * "⠑⢄" * Char(10240) * "\n" * Char(10240)^4 * ""
2329-
@test String(take!(io)) == brailleString
2396+
for transform in (identity, adjoint, transpose)
2397+
_show_with_braille_patterns(convert(IOContext, io), transform(A))
2398+
@test String(take!(io)) == brailleString
2399+
end
23302400

23312401
# Issue #30589
23322402
@test repr("text/plain", sparse([true true])) == "1×2 SparseArrays.SparseMatrixCSC{Bool,$Int} with 2 stored entries:\n 1 1"

0 commit comments

Comments
 (0)