-
Notifications
You must be signed in to change notification settings - Fork 15
Expand file tree
/
Copy pathSparseMatrixCSR.jl
More file actions
212 lines (181 loc) · 4.69 KB
/
SparseMatrixCSR.jl
File metadata and controls
212 lines (181 loc) · 4.69 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
module SparseMatrixCSRTests
using Test
using SparseMatricesCSR
using SparseArrays
using LinearAlgebra
function test_csr(Bi,Tv,Ti)
maxnz=10
maxrows=5
maxcols=6
I = rand(Ti(1):Ti(maxrows),maxnz)
J = rand(Ti(1):Ti(maxcols),maxnz)
V = rand(Tv,maxnz)
#I = Ti[1,2,3]
#J = Ti[1,1,3]
#V = Tv[3,4,6]
CSC = sparse(I,J,V)
if Bi == 1
CSR = sparsecsr(I,J,V)
@test CSR == CSC
@test copy(CSR) == CSC
end
CSR = sparsecsr(Val(Bi),I,J,V)
show(IOContext(stdout, :limit=>true, :displaysize=>(10,10)), CSR)
show(IOContext(stdout, :limit=>false), CSR)
@test CSR == CSC
@test copy(CSR) == CSC
@test eltype(CSR) == Tv
@test isa(CSR,SparseMatrixCSR{Bi,Tv,Ti})
for i=1:size(CSR,1)
for j=1:size(CSR,2)
if (i,j) in zip(I,J)
CSR[i,j] = eltype(V)(i+j)
@test CSR[i,j] ≈ eltype(V)(i+j)
else
try
CSR[i,j] = eltype(V)(i+j)
catch e
@test isa(e,ArgumentError)
end
end
end
end
if Ti == Int64
dense = rand(Tv, maxrows,maxcols)
CSR = sparsecsr(Val(Bi), dense)
@test isa(CSR,SparseMatrixCSR{Bi,Tv,Ti})
@test CSR == dense
if Bi == 1
CSR = sparsecsr(dense)
@test CSR == dense
@test isa(CSR,SparseMatrixCSR{Bi,Tv,Ti})
end
end
CSC = sparse(I,J,V,maxrows,maxcols)
if Bi == 1
CSR = sparsecsr(I,J,V,maxrows,maxcols)
@test CSR == CSC
@test copy(CSR) == CSC
end
CSR = sparsecsr(Val(Bi),I,J,V,maxrows,maxcols)
@test CSR == CSC
@test copy(CSR) == CSC
CSR = sparsecsr(Val(Bi),I,J,V,maxrows,maxcols)
@test CSR == CSC
@test copy(CSR) == CSC
CSR2 = convert(typeof(CSR),CSR)
@test CSR2 === CSR
CSR2 = convert(SparseMatrixCSR{0,Float64,Int32},CSR)
@test CSR2 == CSR
@test copy(CSR2) == CSR
CSR2 = convert(SparseMatrixCSR{0,Float64,Int32},CSC)
@test CSR2 == CSC
@test copy(CSR2) == CSC
CSR2 = convert(SparseMatrixCSR{0,Float64,Int32},collect(CSC))
@test CSR2 == CSC
@test copy(CSR2) == CSC
@test size(CSR) == size(CSC)
@test eltype(CSR) == Tv
@test isa(CSR,SparseMatrixCSR{Bi,Tv,Ti})
@test issparse(CSR)
@test getBi(CSR) == Bi
@test getoffset(CSR) == 1-Bi
@test nnz(CSR) == nnz(CSC)
@test length(SparseArrays.nzvalview(CSR)) == nnz(CSR)
@test nonzeros(CSR) === CSR.nzval
@test colvals(CSR) === CSR.colval
i,j,v = findnz(CSR)
csr = sparsecsr(Val(Bi),i,j,v,maxrows,maxcols)
@test csr == CSR
@test count(v->v>0,CSR) == count(v->v>0,nonzeros(CSR))
x = rand(Tv,maxcols)
y = Vector{Tv}(undef,maxrows)
z = Vector{Tv}(undef,maxrows)
mul!(y,CSR,x)
mul!(z,CSC,x)
@test y ≈ z
@test CSR*x ≈ CSC*x
@test CSR'*y ≈ CSC'*z
mul!(y,CSR,x,1,2)
mul!(z,CSC,x,1,2)
@test y ≈ z
out = LinearAlgebra.fillstored!(CSR,3.33)
@test out === CSR
LinearAlgebra.fillstored!(CSC,3.33)
mul!(y,CSR,x)
mul!(z,CSC,x)
@test y ≈ z
@test_throws ArgumentError fill!(CSR2,3.33)
fill!(CSR2, 0)
@test all(iszero, CSR2)
# test constructors
@test CSR == SparseMatrixCSR(CSC)
@test CSR == SparseMatrixCSR(Matrix(CSC))
_CSR = copy(CSR)
out = LinearAlgebra.rmul!(CSR,-1)
@test out === CSR
@test _CSR ≈ -1*CSR
# Test overlong buffers
let A = sparsecsr([1, 2, 3], [1, 2, 3], [4., 5, 6])
push!(A.nzval, 4.0)
push!(A.rowptr, -1)
@test nnz(A) == length(SparseArrays.nzvalview(A)) == 3
@test SparseArrays.nzvalview(A) == [4., 5, 6]
end
# spzeroscsr tests
csc = spzeros(10, 9)
csr = spzeroscsr(10, 9)
@test csc == csr
csc = spzeros(Tv, 10, 9)
csr = spzeroscsr(Tv, 10, 9)
@test csc == csr
csc = spzeros(Tv, Ti, 10, 9)
csr = spzeroscsr(Tv, Ti, 10, 9)
@test csc == csr
csc = spzeros((10, 9))
csr = spzeroscsr((10, 9))
@test csc == csr
csc = spzeros(Tv, (10, 9))
csr = spzeroscsr(Tv, (10, 9))
@test csc == csr
csc = spzeros(Tv, Ti, (10, 9))
csr = spzeroscsr(Tv, Ti, (10, 9))
@test csc == csr
csc = spzeros(I, J)
csr = spzeroscsr(I, J)
@test csc == csr
csc = spzeros(I, J, maxrows, maxcols)
csr = spzeroscsr(I, J, maxrows, maxcols)
@test csc == csr
csc = spzeros(Tv, I, J, maxrows, maxcols)
csr = spzeroscsr(Tv, I, J, maxrows, maxcols)
@test csc == csr
end
function test_lu(Bi,I,J,V)
CSR=sparsecsr(Val(Bi),I,J,V)
CSC=sparse(I,J,V)
x=rand(3)
@test norm(CSR\x-CSC\x) < 1.0e-14
fact=lu(CSR)
lu!(fact,CSR)
y=similar(x)
ldiv!(y,fact,x)
@test norm(y-CSC\x) < 1.0e-14
end
for Bi in (0,1)
for Tv in (Float32,Float64)
for Ti in (Int32,Int64)
test_csr(Bi,Tv,Ti)
end
end
end
if Base.USE_GPL_LIBS # `lu!` requires `SuiteSparse.UMFPACK`
I = [1,1,2,2,2,3,3]
J = [1,2,1,2,3,2,3]
V = [4.0,1.0,-1.0,4.0,1.0,-1.0,4.0]
test_lu(0,I,J,V)
test_lu(1,I,J,V)
else
@warn "Tests run without GPL libraries."
end
end # module