-
Notifications
You must be signed in to change notification settings - Fork 227
/
conversions.jl
682 lines (601 loc) · 27.3 KB
/
conversions.jl
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
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
export sort_csc, sort_csr, sort_coo
adjtrans_wrappers = ((identity, identity),
(M -> :(Transpose{T, <:$M}), M -> :(_sptranspose(parent($M)))),
(M -> :(Adjoint{T, <:$M}), M -> :(_spadjoint(parent($M)))))
# conversion routines between different sparse and dense storage formats
"""
sparse(x::DenseCuMatrix; fmt=:csc)
sparse(I::CuVector, J::CuVector, V::CuVector, [m, n]; fmt=:csc)
Return a sparse cuda matrix, with type determined by `fmt`.
Possible formats are :csc, :csr, :bsr, and :coo.
"""
function SparseArrays.sparse(x::DenseCuMatrix; fmt=:csc)
if fmt == :csc
return CuSparseMatrixCSC(x)
elseif fmt == :csr
return CuSparseMatrixCSR(x)
elseif fmt == :bsr
return CuSparseMatrixBSR(x)
elseif fmt == :coo
return CuSparseMatrixCOO(x)
else
error("Format :$fmt not available, use :csc, :csr, :bsr or :coo.")
end
end
function SparseArrays.sparse(I::CuVector, J::CuVector, V::CuVector, args...; kwargs...)
sparse(Cint.(I), Cint.(J), V, args...; kwargs...)
end
function SparseArrays.sparse(I::CuVector{Cint}, J::CuVector{Cint}, V::CuVector{Tv},
m=maximum(I), n=maximum(J);
fmt=:csc, combine=nothing) where Tv
# we use COO as an intermediate format, as it's easy to construct from I/J/V.
coo = CuSparseMatrixCOO{Tv}(I, J, V, (m, n))
# find groups of values that correspond to the same position in the matrix.
# if there's no duplicates, `groups` will just be a vector of ones.
# otherwise, it will contain the number of duplicates for each group,
# with subsequent values that are part of the group set to zero.
coo = sort_coo(coo, 'R')
groups = similar(I, Int)
function find_groups(groups, I, J)
i = threadIdx().x + (blockIdx().x - 1) * blockDim().x
if i > length(groups)
return
end
len = 0
# check if we're at the start of a new group
@inbounds if i == 1 || I[i] != I[i-1] || J[i] != J[i-1]
len = 1
while i+len <= length(groups) && I[i] == I[i+len] && J[i] == J[i+len]
len += 1
end
end
@inbounds groups[i] = len
return
end
kernel = @cuda launch=false find_groups(groups, coo.rowInd, coo.colInd)
config = launch_configuration(kernel.fun)
threads = min(length(groups), config.threads)
blocks = cld(length(groups), threads)
kernel(groups, coo.rowInd, coo.colInd; threads, blocks)
# if we got any group of more than one element, we need to combine them.
# this may actually not be required, as some CUSPARSE functions can handle
# duplicate entries, but it's not clear which ones do and which ones don't.
# also, to ensure matrix display is correct, combine values eagerly.
ngroups = mapreduce(!iszero, +, groups)
if ngroups != length(groups)
if combine === nothing
combine = if Tv === Bool
|
else
+
end
end
total_lengths = cumsum(groups) # TODO: add and use `scan!(; exclusive=true)`
I = similar(I, ngroups)
J = similar(J, ngroups)
V = similar(V, ngroups)
# use one thread per value, and if it's at the start of a group,
# combine (if needed) all values and update the output vectors.
function combine_groups(groups, total_lengths, oldI, oldJ, oldV, newI, newJ, newV, combine)
i = threadIdx().x + (blockIdx().x - 1) * blockDim().x
if i > length(groups)
return
end
# check if we're at the start of a group
@inbounds if groups[i] != 0
# get an exclusive offset from the inclusive cumsum
offset = total_lengths[i] - groups[i] + 1
# copy values
newI[i] = oldI[offset]
newJ[i] = oldJ[offset]
newV[i] = if groups[i] == 1
oldV[offset]
else
# combine all values in the group
val = oldV[offset]
j = 1
while j < groups[i]
val = combine(val, oldV[offset+j])
j += 1
end
val
end
end
return
end
kernel = @cuda launch=false combine_groups(groups, total_lengths, coo.rowInd, coo.colInd, coo.nzVal, I, J, V, combine)
config = launch_configuration(kernel.fun)
threads = min(length(groups), config.threads)
blocks = cld(length(groups), threads)
kernel(groups, total_lengths, coo.rowInd, coo.colInd, coo.nzVal, I, J, V, combine; threads, blocks)
coo = CuSparseMatrixCOO{Tv}(I, J, V, (m, n))
end
if fmt == :coo
return coo
elseif fmt == :csc
return CuSparseMatrixCSC(coo)
elseif fmt == :csr
return CuSparseMatrixCSR(coo)
else
error("Format :$fmt not available, use :csc, :csr, or :coo.")
end
end
for (wrapa, unwrapa) in adjtrans_wrappers
for SparseMatrixType in (:(CuSparseMatrixCSC{T}), :(CuSparseMatrixCSR{T}), :(CuSparseMatrixCOO{T}))
TypeA = wrapa(SparseMatrixType)
@eval SparseArrays.sparse(A::$TypeA) where {T} = $(unwrapa(:A))
end
end
function sort_csc(A::CuSparseMatrixCSC{Tv,Ti}, index::SparseChar='O') where {Tv,Ti}
m,n = size(A)
perm = CuArray{Ti}(undef, nnz(A))
cusparseCreateIdentityPermutation(handle(), nnz(A), perm)
descA = CuMatrixDescriptor('G', 'L', 'N', index)
sorted_colPtr = copy(A.colPtr)
sorted_rowVal = copy(A.rowVal)
function bufferSize()
out = Ref{Csize_t}()
cusparseXcscsort_bufferSizeExt(handle(), m, n, nnz(A), A.colPtr, A.rowVal, out)
return out[]
end
with_workspace(bufferSize) do buffer
cusparseXcscsort(handle(), m, n, nnz(A), descA, sorted_colPtr, sorted_rowVal, perm, buffer)
end
perm .+= one(Ti)
sorted_nzVal = A.nzVal[perm]
CUDA.unsafe_free!(perm)
CuSparseMatrixCSC{Tv,Ti}(sorted_colPtr, sorted_rowVal, sorted_nzVal, size(A))
end
function sort_csr(A::CuSparseMatrixCSR{Tv,Ti}, index::SparseChar='O') where {Tv,Ti}
m,n = size(A)
perm = CuArray{Ti}(undef, nnz(A))
cusparseCreateIdentityPermutation(handle(), nnz(A), perm)
descA = CuMatrixDescriptor('G', 'L', 'N', index)
sorted_rowPtr = copy(A.rowPtr)
sorted_colVal = copy(A.colVal)
function bufferSize()
out = Ref{Csize_t}()
cusparseXcsrsort_bufferSizeExt(handle(), m, n, nnz(A), A.rowPtr, A.colVal, out)
return out[]
end
with_workspace(bufferSize) do buffer
cusparseXcsrsort(handle(), m, n, nnz(A), descA, sorted_rowPtr, sorted_colVal, perm, buffer)
end
perm .+= one(Ti)
sorted_nzVal = A.nzVal[perm]
CUDA.unsafe_free!(perm)
CuSparseMatrixCSR{Tv,Ti}(sorted_rowPtr, sorted_colVal, sorted_nzVal, size(A))
end
function sort_coo(A::CuSparseMatrixCOO{Tv,Ti}, type::SparseChar='R') where {Tv,Ti}
type == 'R' || type == 'C' || throw(ArgumentError("type=$type was used and only type='R' and type='C' are supported."))
m,n = size(A)
perm = CuArray{Ti}(undef, nnz(A))
cusparseCreateIdentityPermutation(handle(), nnz(A), perm)
sorted_rowInd = copy(A.rowInd)
sorted_colInd = copy(A.colInd)
function bufferSize()
# It seems that in some cases `out` is not updated
# and we have the following error in the tests:
# "Out of GPU memory trying to allocate 127.781 TiB".
# We set 0 as default value to avoid it.
out = Ref{Csize_t}(0)
cusparseXcoosort_bufferSizeExt(handle(), m, n, nnz(A), A.rowInd, A.colInd, out)
return out[]
end
with_workspace(bufferSize) do buffer
type == 'R' && cusparseXcoosortByRow(handle(), m, n, nnz(A), sorted_rowInd, sorted_colInd, perm, buffer)
type == 'C' && cusparseXcoosortByColumn(handle(), m, n, nnz(A), sorted_rowInd, sorted_colInd, perm, buffer)
end
perm .+= one(Ti)
sorted_nzVal = A.nzVal[perm]
CUDA.unsafe_free!(perm)
CuSparseMatrixCOO{Tv,Ti}(sorted_rowInd, sorted_colInd, sorted_nzVal, size(A))
end
for (bname, fname, pname, elty) in ((:cusparseSpruneCsr2csr_bufferSizeExt, :cusparseSpruneCsr2csrNnz, :cusparseSpruneCsr2csr, :Float32),
(:cusparseDpruneCsr2csr_bufferSizeExt, :cusparseDpruneCsr2csrNnz, :cusparseDpruneCsr2csr, :Float64))
@eval begin
function prune(A::CuSparseMatrixCSR{$elty}, threshold::Number, index::SparseChar)
m, n = size(A)
descA = CuMatrixDescriptor('G', 'L', 'N', index)
descC = CuMatrixDescriptor('G', 'L', 'N', index)
rowPtrC = CuVector{Int32}(undef, m+1)
local colValC, nzValC
function bufferSize()
out = Ref{Csize_t}()
$bname(handle(), m, n, nnz(A), descA, nonzeros(A), A.rowPtr, A.colVal,
Ref{$elty}(threshold), descC, CuPtr{$elty}(CU_NULL), rowPtrC, CuPtr{Int32}(CU_NULL), out)
return out[]
end
with_workspace(bufferSize) do buffer
nnzTotal = Ref{Cint}()
$fname(handle(), m, n, nnz(A), descA, nonzeros(A), A.rowPtr, A.colVal,
Ref{$elty}(threshold), descC, rowPtrC, nnzTotal, buffer)
colValC = CuVector{Int32}(undef, nnzTotal[])
nzValC = CuVector{$elty}(undef, nnzTotal[])
$pname(handle(), m, n, nnz(A), descA, nonzeros(A), A.rowPtr, A.colVal,
Ref{$elty}(threshold), descC, nzValC, rowPtrC, colValC, buffer)
end
return CuSparseMatrixCSR(rowPtrC, colValC, nzValC, (m, n))
end
function prune(A::CuSparseMatrixCSC{$elty}, threshold::Number, index::SparseChar)
m, n = size(A)
descA = CuMatrixDescriptor('G', 'L', 'N', index)
descC = CuMatrixDescriptor('G', 'L', 'N', index)
colPtrC = CuVector{Int32}(undef, n+1)
local rowValC, nzValC
function bufferSize()
out = Ref{Csize_t}()
$bname(handle(), n, m, nnz(A), descA, nonzeros(A), A.colPtr, A.rowVal,
Ref{$elty}(threshold), descC, CuPtr{$elty}(CU_NULL), colPtrC, CuPtr{Int32}(CU_NULL), out)
return out[]
end
with_workspace(bufferSize) do buffer
nnzTotal = Ref{Cint}()
$fname(handle(), n, m, nnz(A), descA, nonzeros(A), A.colPtr, A.rowVal,
Ref{$elty}(threshold), descC, colPtrC, nnzTotal, buffer)
rowValC = CuVector{Int32}(undef, nnzTotal[])
nzValC = CuVector{$elty}(undef, nnzTotal[])
$pname(handle(), n, m, nnz(A), descA, nonzeros(A), A.colPtr, A.rowVal,
Ref{$elty}(threshold), descC, nzValC, colPtrC, rowValC, buffer)
end
return CuSparseMatrixCSC(colPtrC, rowValC, nzValC, (m, n))
end
end
end
## CSR to CSC
function CuSparseMatrixCSC{T}(S::Transpose{T, <:CuSparseMatrixCSR{T}}) where T
csr = parent(S)
return CuSparseMatrixCSC{T}(csr.rowPtr, csr.colVal, csr.nzVal, size(csr))
end
function CuSparseMatrixCSR{T}(S::Transpose{T, <:CuSparseMatrixCSC{T}}) where T
csc = parent(S)
return CuSparseMatrixCSR{T}(csc.colPtr, csc.rowVal, csc.nzVal, size(csc))
end
function CuSparseMatrixCSC{T}(S::Adjoint{T, <:CuSparseMatrixCSR{T}}) where {T <: Real}
csr = parent(S)
return CuSparseMatrixCSC{T}(csr.rowPtr, csr.colVal, csr.nzVal, size(csr))
end
function CuSparseMatrixCSR{T}(S::Adjoint{T, <:CuSparseMatrixCSC{T}}) where {T <: Real}
csc = parent(S)
return CuSparseMatrixCSR{T}(csc.colPtr, csc.rowVal, csc.nzVal, size(csc))
end
function CuSparseMatrixCSC{T}(S::Adjoint{T, <:CuSparseMatrixCSR{T}}) where {T <: Complex}
csr = parent(S)
return CuSparseMatrixCSC{T}(csr.rowPtr, csr.colVal, conj.(csr.nzVal), size(csr))
end
function CuSparseMatrixCSR{T}(S::Adjoint{T, <:CuSparseMatrixCSC{T}}) where {T <: Complex}
csc = parent(S)
return CuSparseMatrixCSR{T}(csc.colPtr, csc.rowVal, conj.(csc.nzVal), size(csc))
end
for SparseMatrixType in (:CuSparseMatrixCSC, :CuSparseMatrixCSR, :CuSparseMatrixCOO)
@eval begin
$SparseMatrixType(S::Diagonal{Tv, <:AbstractVector}) where {Tv} = $SparseMatrixType(cu(S))
$SparseMatrixType(S::Diagonal{Tv, <:CuArray}) where Tv = $SparseMatrixType{Tv}(S)
$SparseMatrixType{Tv}(S::Diagonal) where {Tv} = $SparseMatrixType{Tv, Cint}(S)
end
if SparseMatrixType == :CuSparseMatrixCOO
@eval function $SparseMatrixType{Tv, Ti}(S::Diagonal) where {Tv, Ti}
m = size(S, 1)
return $SparseMatrixType{Tv, Ti}(CuVector(1:m), CuVector(1:m), convert(CuVector{Tv}, S.diag), (m, m))
end
else
@eval function $SparseMatrixType{Tv, Ti}(S::Diagonal) where {Tv, Ti}
m = size(S, 1)
return $SparseMatrixType{Tv, Ti}(CuVector(1:(m+1)), CuVector(1:m), convert(CuVector{Tv}, S.diag), (m, m))
end
end
end
# by flipping rows and columns, we can use that to get CSC to CSR too
for elty in (:Float32, :Float64, :ComplexF32, :ComplexF64)
@eval begin
function CuSparseMatrixCSC{$elty}(csr::CuSparseMatrixCSR{$elty}; index::SparseChar='O', action::cusparseAction_t=CUSPARSE_ACTION_NUMERIC, algo::cusparseCsr2CscAlg_t=CUSPARSE_CSR2CSC_ALG1)
m,n = size(csr)
colPtr = CUDA.zeros(Cint, n+1)
rowVal = CUDA.zeros(Cint, nnz(csr))
nzVal = CUDA.zeros($elty, nnz(csr))
function bufferSize()
out = Ref{Csize_t}(1)
cusparseCsr2cscEx2_bufferSize(handle(), m, n, nnz(csr), nonzeros(csr),
csr.rowPtr, csr.colVal, nzVal, colPtr, rowVal,
$elty, action, index, algo, out)
return out[]
end
with_workspace(bufferSize) do buffer
cusparseCsr2cscEx2(handle(), m, n, nnz(csr), nonzeros(csr),
csr.rowPtr, csr.colVal, nzVal, colPtr, rowVal,
$elty, action, index, algo, buffer)
end
CuSparseMatrixCSC(colPtr,rowVal,nzVal,size(csr))
end
function CuSparseMatrixCSR{$elty}(csc::CuSparseMatrixCSC{$elty}; index::SparseChar='O', action::cusparseAction_t=CUSPARSE_ACTION_NUMERIC, algo::cusparseCsr2CscAlg_t=CUSPARSE_CSR2CSC_ALG1)
m,n = size(csc)
rowPtr = CUDA.zeros(Cint,m+1)
colVal = CUDA.zeros(Cint,nnz(csc))
nzVal = CUDA.zeros($elty,nnz(csc))
function bufferSize()
out = Ref{Csize_t}(1)
cusparseCsr2cscEx2_bufferSize(handle(), n, m, nnz(csc), nonzeros(csc),
csc.colPtr, rowvals(csc), nzVal, rowPtr, colVal,
$elty, action, index, algo, out)
return out[]
end
with_workspace(bufferSize) do buffer
cusparseCsr2cscEx2(handle(), n, m, nnz(csc), nonzeros(csc),
csc.colPtr, rowvals(csc), nzVal, rowPtr, colVal,
$elty, action, index, algo, buffer)
end
CuSparseMatrixCSR(rowPtr,colVal,nzVal,size(csc))
end
end
end
# implement Float16 conversions using wider types
# TODO: Float16 is sometimes natively supported
for (elty, welty) in ((:Float16, :Float32),
(:ComplexF16, :ComplexF32))
@eval begin
function CuSparseMatrixCSC{$elty}(csr::CuSparseMatrixCSR{$elty}; index::SparseChar='O', action::cusparseAction_t=CUSPARSE_ACTION_NUMERIC, algo::cusparseCsr2CscAlg_t=CUSPARSE_CSR2CSC_ALG1)
m,n = size(csr)
colPtr = CUDA.zeros(Cint, n+1)
rowVal = CUDA.zeros(Cint, nnz(csr))
nzVal = CUDA.zeros($elty, nnz(csr))
if $elty == Float16 #broken for ComplexF16?
function bufferSize()
out = Ref{Csize_t}(1)
cusparseCsr2cscEx2_bufferSize(handle(), m, n, nnz(csr), nonzeros(csr),
csr.rowPtr, csr.colVal, nzVal, colPtr, rowVal,
$elty, action, index, algo, out)
return out[]
end
with_workspace(bufferSize) do buffer
cusparseCsr2cscEx2(handle(), m, n, nnz(csr), nonzeros(csr),
csr.rowPtr, csr.colVal, nzVal, colPtr, rowVal,
$elty, action, index, algo, buffer)
end
return CuSparseMatrixCSC(colPtr,rowVal,nzVal,size(csr))
else
wide_csr = CuSparseMatrixCSR(csr.rowPtr, csr.colVal, convert(CuVector{$welty}, nonzeros(csr)), size(csr))
wide_csc = CuSparseMatrixCSC(wide_csr)
return CuSparseMatrixCSC(wide_csc.colPtr, wide_csc.rowVal, convert(CuVector{$elty}, nonzeros(wide_csc)), size(wide_csc))
end
end
function CuSparseMatrixCSR{$elty}(csc::CuSparseMatrixCSC{$elty}; index::SparseChar='O', action::cusparseAction_t=CUSPARSE_ACTION_NUMERIC, algo::cusparseCsr2CscAlg_t=CUSPARSE_CSR2CSC_ALG1)
m,n = size(csc)
rowPtr = CUDA.zeros(Cint,m+1)
colVal = CUDA.zeros(Cint,nnz(csc))
nzVal = CUDA.zeros($elty,nnz(csc))
if $elty == Float16 #broken for ComplexF16?
function bufferSize()
out = Ref{Csize_t}(1)
cusparseCsr2cscEx2_bufferSize(handle(), n, m, nnz(csc), nonzeros(csc),
csc.colPtr, rowvals(csc), nzVal, rowPtr, colVal,
$elty, action, index, algo, out)
return out[]
end
with_workspace(bufferSize) do buffer
cusparseCsr2cscEx2(handle(), n, m, nnz(csc), nonzeros(csc),
csc.colPtr, rowvals(csc), nzVal, rowPtr, colVal,
$elty, action, index, algo, buffer)
end
return CuSparseMatrixCSR(rowPtr,colVal,nzVal,size(csc))
else
wide_csc = CuSparseMatrixCSC(csc.colPtr, csc.rowVal, convert(CuVector{$welty}, nonzeros(csc)), size(csc))
wide_csr = CuSparseMatrixCSR(wide_csc)
return CuSparseMatrixCSR(wide_csr.rowPtr, wide_csr.colVal, convert(CuVector{$elty}, nonzeros(wide_csr)), size(wide_csr))
end
end
end
end
# implement Int conversions using reinterpreted Float
for (elty, felty) in ((:Int16, :Float16),
(:Int32, :Float32),
(:Int64, :Float64),
(:Int128, :ComplexF64))
@eval begin
function CuSparseMatrixCSR{$elty}(csc::CuSparseMatrixCSC{$elty})
csc_compat = CuSparseMatrixCSC(
csc.colPtr,
csc.rowVal,
reinterpret($felty, csc.nzVal),
size(csc)
)
csr_compat = CuSparseMatrixCSR(csc_compat)
CuSparseMatrixCSR(
csr_compat.rowPtr,
csr_compat.colVal,
reinterpret($elty, csr_compat.nzVal),
size(csr_compat)
)
end
function CuSparseMatrixCSC{$elty}(csr::CuSparseMatrixCSR{$elty})
csr_compat = CuSparseMatrixCSR(
csr.rowPtr,
csr.colVal,
reinterpret($felty, csr.nzVal),
size(csr)
)
csc_compat = CuSparseMatrixCSC(csr_compat)
CuSparseMatrixCSC(
csc_compat.colPtr,
csc_compat.rowVal,
reinterpret($elty, csc_compat.nzVal),
size(csc_compat)
)
end
end
end
## CuSparseVector to CuVector
CuVector(x::CuSparseVector{T}) where {T} = CuVector{T}(x)
function CuVector{T}(sv::CuSparseVector{T}) where {T}
n = length(sv)
dv = CUDA.zeros(T, n)
scatter!(dv, sv, 'O')
end
## CSR to BSR and vice-versa
for (fname,elty) in ((:cusparseScsr2bsr, :Float32),
(:cusparseDcsr2bsr, :Float64),
(:cusparseCcsr2bsr, :ComplexF32),
(:cusparseZcsr2bsr, :ComplexF64))
@eval begin
function CuSparseMatrixBSR{$elty}(csr::CuSparseMatrixCSR{$elty}, blockDim::Integer;
dir::SparseChar='R', index::SparseChar='O',
indc::SparseChar='O')
m,n = size(csr)
nnz_ref = Ref{Cint}(1)
mb = cld(m, blockDim)
nb = cld(n, blockDim)
bsrRowPtr = CUDA.zeros(Cint,mb + 1)
cudesca = CuMatrixDescriptor('G', 'L', 'N', index)
cudescc = CuMatrixDescriptor('G', 'L', 'N', indc)
cusparseXcsr2bsrNnz(handle(), dir, m, n, cudesca, csr.rowPtr,
csr.colVal, blockDim, cudescc, bsrRowPtr, nnz_ref)
(nnz_ref[] > mb * nb) && error("The number of nonzero blocks of the BSR matrix is incorrect.")
bsrNzVal = CUDA.zeros($elty, nnz_ref[] * blockDim * blockDim )
bsrColInd = CUDA.zeros(Cint, nnz_ref[])
$fname(handle(), dir, m, n,
cudesca, nonzeros(csr), csr.rowPtr, csr.colVal,
blockDim, cudescc, bsrNzVal, bsrRowPtr,
bsrColInd)
CuSparseMatrixBSR{$elty}(bsrRowPtr, bsrColInd, bsrNzVal, size(csr), blockDim, dir, nnz_ref[])
end
end
end
for (fname,elty) in ((:cusparseSbsr2csr, :Float32),
(:cusparseDbsr2csr, :Float64),
(:cusparseCbsr2csr, :ComplexF32),
(:cusparseZbsr2csr, :ComplexF64))
@eval begin
function CuSparseMatrixCSR{$elty}(bsr::CuSparseMatrixBSR{$elty};
index::SparseChar='O', indc::SparseChar='O')
m,n = size(bsr)
mb = cld(m, bsr.blockDim)
nb = cld(n, bsr.blockDim)
cudesca = CuMatrixDescriptor('G', 'L', 'N', index)
cudescc = CuMatrixDescriptor('G', 'L', 'N', indc)
csrRowPtr = CUDA.zeros(Cint, m + 1)
csrColInd = CUDA.zeros(Cint, nnz(bsr))
csrNzVal = CUDA.zeros($elty, nnz(bsr))
$fname(handle(), bsr.dir, mb, nb,
cudesca, nonzeros(bsr), bsr.rowPtr, bsr.colVal,
bsr.blockDim, cudescc, csrNzVal, csrRowPtr,
csrColInd)
# XXX: the size here may not match the expected size, when the matrix dimension
# is not a multiple of the block dimension!
CuSparseMatrixCSR(csrRowPtr, csrColInd, csrNzVal, (mb*bsr.blockDim, nb*bsr.blockDim))
end
end
end
# implement Int conversions using reinterpreted Float
for (elty, felty) in ((:Int16, :Float16),
(:Int32, :Float32),
(:Int64, :Float64),
(:Int128, :ComplexF64))
@eval begin
function CuSparseMatrixCSR{$elty}(bsr::CuSparseMatrixBSR{$elty})
bsr_compat = CuSparseMatrixBSR(
bsr.rowPtr,
bsr.colVal,
reinterpret($felty, bsr.nzVal),
bsr.blockDim,
bsr.dir,
bsr.nnzb,
size(bsr)
)
csr_compat = CuSparseMatrixCSR(bsr_compat)
CuSparseMatrixCSR(
csr_compat.rowPtr,
csr_compat.colVal,
reinterpret($elty, csr_compat.nzVal),
size(csr_compat)
)
end
function CuSparseMatrixBSR{$elty}(csr::CuSparseMatrixCSR{$elty}, blockDim)
csr_compat = CuSparseMatrixCSR(
csr.rowPtr,
csr.colVal,
reinterpret($felty, csr.nzVal),
size(csr)
)
bsr_compat = CuSparseMatrixBSR(csr_compat, blockDim)
CuSparseMatrixBSR(
bsr_compat.rowPtr,
bsr_compat.colVal,
reinterpret($elty, bsr_compat.nzVal),
bsr_compat.blockDim,
bsr_compat.dir,
bsr_compat.nnzb,
size(bsr_compat)
)
end
end
end
## CSR to COO and vice-versa
function CuSparseMatrixCSR{Tv}(coo::CuSparseMatrixCOO{Tv}; index::SparseChar='O') where {Tv}
m,n = size(coo)
coo = sort_coo(coo, 'R')
csrRowPtr = CuVector{Cint}(undef, m+1)
cusparseXcoo2csr(handle(), coo.rowInd, nnz(coo), m, csrRowPtr, index)
CuSparseMatrixCSR{Tv}(csrRowPtr, coo.colInd, nonzeros(coo), size(coo))
end
function CuSparseMatrixCOO{Tv}(csr::CuSparseMatrixCSR{Tv}; index::SparseChar='O') where {Tv}
m,n = size(csr)
cooRowInd = CuVector{Cint}(undef, nnz(csr))
cusparseXcsr2coo(handle(), csr.rowPtr, nnz(csr), m, cooRowInd, index)
CuSparseMatrixCOO{Tv}(cooRowInd, csr.colVal, nonzeros(csr), size(csr))
end
### CSC to COO and viceversa
function CuSparseMatrixCSC{Tv}(coo::CuSparseMatrixCOO{Tv}; index::SparseChar='O') where {Tv}
m,n = size(coo)
coo = sort_coo(coo, 'C')
cscColPtr = CuVector{Cint}(undef, n+1)
cusparseXcoo2csr(handle(), coo.colInd, nnz(coo), n, cscColPtr, index)
CuSparseMatrixCSC{Tv}(cscColPtr, coo.rowInd, nonzeros(coo), size(coo))
end
function CuSparseMatrixCOO{Tv}(csc::CuSparseMatrixCSC{Tv}; index::SparseChar='O') where {Tv}
m,n = size(csc)
cooColInd = CuVector{Cint}(undef, nnz(csc))
cusparseXcsr2coo(handle(), csc.colPtr, nnz(csc), n, cooColInd, index)
coo = CuSparseMatrixCOO{Tv}(csc.rowVal, cooColInd, nonzeros(csc), size(csc))
coo = sort_coo(coo, 'R')
end
### BSR to COO and viceversa
CuSparseMatrixBSR(coo::CuSparseMatrixCOO, blockdim) = CuSparseMatrixBSR(CuSparseMatrixCSR(coo), blockdim) # no direct conversion
CuSparseMatrixCOO(bsr::CuSparseMatrixBSR) = CuSparseMatrixCOO(CuSparseMatrixCSR(bsr)) # no direct conversion
### BSR to CSC and viceversa
CuSparseMatrixBSR(csc::CuSparseMatrixCSC, blockdim) = CuSparseMatrixBSR(CuSparseMatrixCSR(csc), blockdim) # no direct conversion
CuSparseMatrixCSC(bsr::CuSparseMatrixBSR) = CuSparseMatrixCSC(CuSparseMatrixCSR(bsr)) # no direct conversion
## sparse to dense, and vice-versa
function CUDA.CuMatrix{T}(csr::CuSparseMatrixCSR{T}; index::SparseChar='O') where {T}
denseA = sparsetodense(csr, index)
return denseA
end
function CUDA.CuMatrix{T}(csc::CuSparseMatrixCSC{T}; index::SparseChar='O') where {T}
denseA = sparsetodense(csc, index)
return denseA
end
Base.copyto!(dest::Matrix{T}, src::AbstractCuSparseMatrix{T}) where T = copyto!(dest, CuMatrix{T}(src))
function CuSparseMatrixCSR(A::CuMatrix{T}; index::SparseChar='O', sorted::Bool=false) where {T}
csr = densetosparse(A, :csr, index)
csr = sorted ? sort_csr(csr, index) : csr
return csr
end
function CuSparseMatrixCSC(A::CuMatrix{T}; index::SparseChar='O', sorted::Bool=false) where {T}
csc = densetosparse(A, :csc, index)
csc = sorted ? sort_csc(csc, index) : csc
return csc
end
function CUDA.CuMatrix{T}(bsr::CuSparseMatrixBSR{T}; index::SparseChar='O',
indc::SparseChar='O') where {T}
CuMatrix{T}(CuSparseMatrixCSR{T}(bsr; index, indc))
end
function CuSparseMatrixBSR(A::CuMatrix, blockDim::Integer=gcd(size(A)...); index::SparseChar='O')
m,n = size(A)
# csr.colVal should be sorted if we want to use "csr2bsr" routines.
csr = CuSparseMatrixCSR(A; index, sorted=true)
CuSparseMatrixBSR(csr, blockDim)
end
function CUDA.CuMatrix{T}(coo::CuSparseMatrixCOO{T}; index::SparseChar='O') where {T}
sparsetodense(coo, index)
end
function CuSparseMatrixCOO(A::CuMatrix{T}; index::SparseChar='O') where {T}
densetosparse(A, :coo, index)
end