-
-
Notifications
You must be signed in to change notification settings - Fork 5.5k
/
sparsematrix.jl
3559 lines (3142 loc) · 122 KB
/
sparsematrix.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
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# This file is a part of Julia. License is MIT: https://julialang.org/license
# Compressed sparse columns data structure
# Assumes that no zeros are stored in the data structure
# Assumes that row values in rowval for each column are sorted
# issorted(rowval[colptr[i]:(colptr[i+1]-1)]) == true
"""
SparseMatrixCSC{Tv,Ti<:Integer} <: AbstractSparseMatrix{Tv,Ti}
Matrix type for storing sparse matrices in the
[Compressed Sparse Column](@ref man-csc) format.
"""
struct SparseMatrixCSC{Tv,Ti<:Integer} <: AbstractSparseMatrix{Tv,Ti}
m::Int # Number of rows
n::Int # Number of columns
colptr::Vector{Ti} # Column i is in colptr[i]:(colptr[i+1]-1)
rowval::Vector{Ti} # Row indices of stored values
nzval::Vector{Tv} # Stored values, typically nonzeros
function SparseMatrixCSC{Tv,Ti}(m::Integer, n::Integer, colptr::Vector{Ti}, rowval::Vector{Ti},
nzval::Vector{Tv}) where {Tv,Ti<:Integer}
@noinline throwsz(str, lbl, k) =
throw(ArgumentError("number of $str ($lbl) must be ≥ 0, got $k"))
m < 0 && throwsz("rows", 'm', m)
n < 0 && throwsz("columns", 'n', n)
new(Int(m), Int(n), colptr, rowval, nzval)
end
end
function SparseMatrixCSC(m::Integer, n::Integer, colptr::Vector, rowval::Vector, nzval::Vector)
Tv = eltype(nzval)
Ti = promote_type(eltype(colptr), eltype(rowval))
SparseMatrixCSC{Tv,Ti}(m, n, colptr, rowval, nzval)
end
size(S::SparseMatrixCSC) = (S.m, S.n)
# Define an alias for views of a SparseMatrixCSC which include all rows and a unit range of the columns.
# Also define a union of SparseMatrixCSC and this view since many methods can be defined efficiently for
# this union by extracting the fields via the get function: getcolptr, getrowval, and getnzval. The key
# insight is that getcolptr on a SparseMatrixCSCView returns an offset view of the colptr of the
# underlying SparseMatrixCSC
const SparseMatrixCSCView{Tv,Ti} =
SubArray{Tv,2,SparseMatrixCSC{Tv,Ti},
Tuple{Base.Slice{Base.OneTo{Int}},I}} where {I<:AbstractUnitRange}
const SparseMatrixCSCUnion{Tv,Ti} = Union{SparseMatrixCSC{Tv,Ti}, SparseMatrixCSCView{Tv,Ti}}
getcolptr(S::SparseMatrixCSC) = S.colptr
getcolptr(S::SparseMatrixCSCView) = view(S.parent.colptr, first(axes(S, 2)):(last(axes(S, 2)) + 1))
getrowval(S::SparseMatrixCSC) = S.rowval
getrowval(S::SparseMatrixCSCView) = S.parent.rowval
getnzval( S::SparseMatrixCSC) = S.nzval
getnzval( S::SparseMatrixCSCView) = S.parent.nzval
nzvalview(S::SparseMatrixCSC) = view(S.nzval, 1:nnz(S))
"""
nnz(A)
Returns the number of stored (filled) elements in a sparse array.
# Examples
```jldoctest
julia> A = sparse(2I, 3, 3)
3×3 SparseMatrixCSC{Int64,Int64} with 3 stored entries:
[1, 1] = 2
[2, 2] = 2
[3, 3] = 2
julia> nnz(A)
3
```
"""
nnz(S::SparseMatrixCSC) = Int(S.colptr[S.n + 1] - 1)
nnz(S::ReshapedArray{T,1,<:SparseMatrixCSC}) where T = nnz(parent(S))
count(pred, S::SparseMatrixCSC) = count(pred, nzvalview(S)) + pred(zero(eltype(S)))*(prod(size(S)) - nnz(S))
"""
nonzeros(A)
Return a vector of the structural nonzero values in sparse array `A`. This
includes zeros that are explicitly stored in the sparse array. The returned
vector points directly to the internal nonzero storage of `A`, and any
modifications to the returned vector will mutate `A` as well. See
[`rowvals`](@ref) and [`nzrange`](@ref).
# Examples
```jldoctest
julia> A = sparse(2I, 3, 3)
3×3 SparseMatrixCSC{Int64,Int64} with 3 stored entries:
[1, 1] = 2
[2, 2] = 2
[3, 3] = 2
julia> nonzeros(A)
3-element Array{Int64,1}:
2
2
2
```
"""
nonzeros(S::SparseMatrixCSC) = S.nzval
"""
rowvals(A::SparseMatrixCSC)
Return a vector of the row indices of `A`. Any modifications to the returned
vector will mutate `A` as well. Providing access to how the row indices are
stored internally can be useful in conjunction with iterating over structural
nonzero values. See also [`nonzeros`](@ref) and [`nzrange`](@ref).
# Examples
```jldoctest
julia> A = sparse(2I, 3, 3)
3×3 SparseMatrixCSC{Int64,Int64} with 3 stored entries:
[1, 1] = 2
[2, 2] = 2
[3, 3] = 2
julia> rowvals(A)
3-element Array{Int64,1}:
1
2
3
```
"""
rowvals(S::SparseMatrixCSC) = S.rowval
"""
nzrange(A::SparseMatrixCSC, col::Integer)
Return the range of indices to the structural nonzero values of a sparse matrix
column. In conjunction with [`nonzeros`](@ref) and
[`rowvals`](@ref), this allows for convenient iterating over a sparse matrix :
A = sparse(I,J,V)
rows = rowvals(A)
vals = nonzeros(A)
m, n = size(A)
for i = 1:n
for j in nzrange(A, i)
row = rows[j]
val = vals[j]
# perform sparse wizardry...
end
end
"""
nzrange(S::SparseMatrixCSC, col::Integer) = S.colptr[col]:(S.colptr[col+1]-1)
function Base.show(io::IO, ::MIME"text/plain", S::SparseMatrixCSC)
xnnz = nnz(S)
print(io, S.m, "×", S.n, " ", typeof(S), " with ", xnnz, " stored ",
xnnz == 1 ? "entry" : "entries")
if xnnz != 0
print(io, ":")
show(io, S)
end
end
Base.show(io::IO, S::SparseMatrixCSC) = Base.show(convert(IOContext, io), S::SparseMatrixCSC)
function Base.show(io::IOContext, S::SparseMatrixCSC)
if nnz(S) == 0
return show(io, MIME("text/plain"), S)
end
limit::Bool = get(io, :limit, false)
rows = displaysize(io)[1] - 4 # -4 from [Prompt, header, newline after elements, new prompt]
will_fit = !limit || rows >= nnz(S) # Will the whole matrix fit when printed?
if rows <= 2 && !will_fit
print(io, "\n \u22ee")
return
end
iob = IOBuffer()
ioc = IOContext(iob, :compact => true)
function _format_line(r, col)
pad = ndigits(max(S.m, S.n))
print(ioc, " [", rpad(S.rowval[r], pad), ", ", lpad(col, pad), "] = ")
if isassigned(S.nzval, Int(r))
show(ioc, S.nzval[r])
else
print(ioc, Base.undef_ref_str)
end
return String(take!(iob))
end
if will_fit
print_count = nnz(S)
else
print_count = div(rows-1, 2)
end
count = 0
for col = 1:S.n, r = nzrange(S, col)
count += 1
print(io, "\n", _format_line(r, col))
count == print_count && break
end
if !will_fit
print(io, "\n \u22ee")
# find the column to start printing in for the last print_count elements
nextcol = searchsortedfirst(S.colptr, nnz(S) - print_count + 1)
for r = (nnz(S) - print_count + 1) : (S.colptr[nextcol] - 1)
print(io, "\n", _format_line(r, nextcol - 1))
end
# print all of the remaining columns
for col = nextcol:S.n, r = nzrange(S, col)
print(io, "\n", _format_line(r, col))
end
end
end
## Reshape
function sparse_compute_reshaped_colptr_and_rowval(colptrS::Vector{Ti}, rowvalS::Vector{Ti},
mS::Int, nS::Int, colptrA::Vector{Ti},
rowvalA::Vector{Ti}, mA::Int, nA::Int) where Ti
lrowvalA = length(rowvalA)
maxrowvalA = (lrowvalA > 0) ? maximum(rowvalA) : zero(Ti)
((length(colptrA) == (nA+1)) && (maximum(colptrA) <= (lrowvalA+1)) && (maxrowvalA <= mA)) || throw(BoundsError())
colptrS[1] = 1
colA = 1
colS = 1
ptr = 1
@inbounds while colA <= nA
offsetA = (colA - 1) * mA
while ptr <= colptrA[colA+1]-1
rowA = rowvalA[ptr]
i = offsetA + rowA - 1
colSn = div(i, mS) + 1
rowS = mod(i, mS) + 1
while colS < colSn
colptrS[colS+1] = ptr
colS += 1
end
rowvalS[ptr] = rowS
ptr += 1
end
colA += 1
end
@inbounds while colS <= nS
colptrS[colS+1] = ptr
colS += 1
end
end
function copy(ra::ReshapedArray{<:Any,2,<:SparseMatrixCSC})
mS,nS = size(ra)
a = parent(ra)
mA,nA = size(a)
numnz = nnz(a)
colptr = similar(a.colptr, nS+1)
rowval = similar(a.rowval)
nzval = copy(a.nzval)
sparse_compute_reshaped_colptr_and_rowval(colptr, rowval, mS, nS, a.colptr, a.rowval, mA, nA)
return SparseMatrixCSC(mS, nS, colptr, rowval, nzval)
end
## Alias detection and prevention
using Base: dataids, unaliascopy
Base.dataids(S::SparseMatrixCSC) = (dataids(S.colptr)..., dataids(S.rowval)..., dataids(S.nzval)...)
Base.unaliascopy(S::SparseMatrixCSC) = typeof(S)(S.m, S.n, unaliascopy(S.colptr), unaliascopy(S.rowval), unaliascopy(S.nzval))
## Constructors
copy(S::SparseMatrixCSC) =
SparseMatrixCSC(S.m, S.n, copy(S.colptr), copy(S.rowval), copy(S.nzval))
function copyto!(A::SparseMatrixCSC, B::SparseMatrixCSC)
# If the two matrices have the same length then all the
# elements in A will be overwritten.
if length(A) == length(B)
resize!(A.nzval, length(B.nzval))
resize!(A.rowval, length(B.rowval))
if size(A) == size(B)
# Simple case: we can simply copy the internal fields of B to A.
copyto!(A.colptr, B.colptr)
copyto!(A.rowval, B.rowval)
else
# This is like a "reshape B into A".
sparse_compute_reshaped_colptr_and_rowval(A.colptr, A.rowval, A.m, A.n, B.colptr, B.rowval, B.m, B.n)
end
else
length(A) >= length(B) || throw(BoundsError())
lB = length(B)
nnzA = nnz(A)
nnzB = nnz(B)
# Up to which col, row, and ptr in rowval/nzval will A be overwritten?
lastmodcolA = div(lB - 1, A.m) + 1
lastmodrowA = mod(lB - 1, A.m) + 1
lastmodptrA = A.colptr[lastmodcolA]
while lastmodptrA < A.colptr[lastmodcolA+1] && A.rowval[lastmodptrA] <= lastmodrowA
lastmodptrA += 1
end
lastmodptrA -= 1
if lastmodptrA >= nnzB
# A will have fewer non-zero elements; unmodified elements are kept at the end.
deleteat!(A.rowval, nnzB+1:lastmodptrA)
deleteat!(A.nzval, nnzB+1:lastmodptrA)
else
# A will have more non-zero elements; unmodified elements are kept at the end.
resize!(A.rowval, nnzB + nnzA - lastmodptrA)
resize!(A.nzval, nnzB + nnzA - lastmodptrA)
copyto!(A.rowval, nnzB+1, A.rowval, lastmodptrA+1, nnzA-lastmodptrA)
copyto!(A.nzval, nnzB+1, A.nzval, lastmodptrA+1, nnzA-lastmodptrA)
end
# Adjust colptr accordingly.
@inbounds for i in 2:length(A.colptr)
A.colptr[i] += nnzB - lastmodptrA
end
sparse_compute_reshaped_colptr_and_rowval(A.colptr, A.rowval, A.m, lastmodcolA-1, B.colptr, B.rowval, B.m, B.n)
end
copyto!(A.nzval, B.nzval)
return A
end
## similar
#
# parent method for similar that preserves stored-entry structure (for when new and old dims match)
function _sparsesimilar(S::SparseMatrixCSC, ::Type{TvNew}, ::Type{TiNew}) where {TvNew,TiNew}
newcolptr = copyto!(similar(S.colptr, TiNew), S.colptr)
newrowval = copyto!(similar(S.rowval, TiNew), S.rowval)
return SparseMatrixCSC(S.m, S.n, newcolptr, newrowval, similar(S.nzval, TvNew))
end
# parent methods for similar that preserves only storage space (for when new and old dims differ)
_sparsesimilar(S::SparseMatrixCSC, ::Type{TvNew}, ::Type{TiNew}, dims::Dims{2}) where {TvNew,TiNew} =
SparseMatrixCSC(dims..., fill(one(TiNew), last(dims)+1), similar(S.rowval, TiNew), similar(S.nzval, TvNew))
# parent method for similar that allocates an empty sparse vector (when new dims are single)
_sparsesimilar(S::SparseMatrixCSC, ::Type{TvNew}, ::Type{TiNew}, dims::Dims{1}) where {TvNew,TiNew} =
SparseVector(dims..., similar(S.rowval, TiNew, 0), similar(S.nzval, TvNew, 0))
#
# The following methods hook into the AbstractArray similar hierarchy. The first method
# covers similar(A[, Tv]) calls, which preserve stored-entry structure, and the latter
# methods cover similar(A[, Tv], shape...) calls, which preserve storage space when the shape
# calls for a two-dimensional result.
similar(S::SparseMatrixCSC{<:Any,Ti}, ::Type{TvNew}) where {Ti,TvNew} = _sparsesimilar(S, TvNew, Ti)
similar(S::SparseMatrixCSC{<:Any,Ti}, ::Type{TvNew}, dims::Union{Dims{1},Dims{2}}) where {Ti,TvNew} =
_sparsesimilar(S, TvNew, Ti, dims)
# The following methods cover similar(A, Tv, Ti[, shape...]) calls, which specify the
# result's index type in addition to its entry type, and aren't covered by the hooks above.
# The calls without shape again preserve stored-entry structure, whereas those with shape
# preserve storage space when the shape calls for a two-dimensional result.
similar(S::SparseMatrixCSC, ::Type{TvNew}, ::Type{TiNew}) where{TvNew,TiNew} =
_sparsesimilar(S, TvNew, TiNew)
similar(S::SparseMatrixCSC, ::Type{TvNew}, ::Type{TiNew}, dims::Union{Dims{1},Dims{2}}) where {TvNew,TiNew} =
_sparsesimilar(S, TvNew, TiNew, dims)
similar(S::SparseMatrixCSC, ::Type{TvNew}, ::Type{TiNew}, m::Integer) where {TvNew,TiNew} =
_sparsesimilar(S, TvNew, TiNew, (m,))
similar(S::SparseMatrixCSC, ::Type{TvNew}, ::Type{TiNew}, m::Integer, n::Integer) where {TvNew,TiNew} =
_sparsesimilar(S, TvNew, TiNew, (m, n))
# converting between SparseMatrixCSC types
SparseMatrixCSC(S::SparseMatrixCSC) = copy(S)
AbstractMatrix{Tv}(A::SparseMatrixCSC) where {Tv} = SparseMatrixCSC{Tv}(A)
SparseMatrixCSC{Tv}(S::SparseMatrixCSC{Tv}) where {Tv} = copy(S)
SparseMatrixCSC{Tv}(S::SparseMatrixCSC) where {Tv} = SparseMatrixCSC{Tv,eltype(S.colptr)}(S)
SparseMatrixCSC{Tv,Ti}(S::SparseMatrixCSC{Tv,Ti}) where {Tv,Ti} = copy(S)
function SparseMatrixCSC{Tv,Ti}(S::SparseMatrixCSC) where {Tv,Ti}
eltypeTicolptr = convert(Vector{Ti}, S.colptr)
eltypeTirowval = convert(Vector{Ti}, S.rowval)
eltypeTvnzval = convert(Vector{Tv}, S.nzval)
return SparseMatrixCSC(S.m, S.n, eltypeTicolptr, eltypeTirowval, eltypeTvnzval)
end
# converting from other matrix types to SparseMatrixCSC (also see sparse())
SparseMatrixCSC(M::Matrix) = sparse(M)
SparseMatrixCSC(M::AbstractMatrix{Tv}) where {Tv} = SparseMatrixCSC{Tv,Int}(M)
SparseMatrixCSC{Tv}(M::AbstractMatrix{Tv}) where {Tv} = SparseMatrixCSC{Tv,Int}(M)
function SparseMatrixCSC{Tv,Ti}(M::AbstractMatrix) where {Tv,Ti}
@assert !has_offset_axes(M)
I = findall(x -> x != 0, M)
eltypeTiI = Ti[i[1] for i in I]
eltypeTiJ = Ti[i[2] for i in I]
eltypeTvV = Tv[M[i] for i in I]
return sparse_IJ_sorted!(eltypeTiI, eltypeTiJ, eltypeTvV, size(M)...)
end
function SparseMatrixCSC{Tv,Ti}(M::StridedMatrix) where {Tv,Ti}
nz = count(t -> t != 0, M)
colptr = zeros(Ti, size(M, 2) + 1)
nzval = Vector{Tv}(undef, nz)
rowval = Vector{Ti}(undef, nz)
colptr[1] = 1
cnt = 1
@inbounds for j in 1:size(M, 2)
for i in 1:size(M, 1)
v = M[i, j]
if v != 0
rowval[cnt] = i
nzval[cnt] = v
cnt += 1
end
end
colptr[j+1] = cnt
end
return SparseMatrixCSC(size(M, 1), size(M, 2), colptr, rowval, nzval)
end
SparseMatrixCSC(M::Adjoint{<:Any,<:SparseMatrixCSC}) = copy(M)
SparseMatrixCSC(M::Transpose{<:Any,<:SparseMatrixCSC}) = copy(M)
SparseMatrixCSC{Tv}(M::Adjoint{Tv,SparseMatrixCSC{Tv}}) where {Tv} = copy(M)
SparseMatrixCSC{Tv}(M::Transpose{Tv,SparseMatrixCSC{Tv}}) where {Tv} = copy(M)
SparseMatrixCSC{Tv,Ti}(M::Adjoint{Tv,SparseMatrixCSC{Tv,Ti}}) where {Tv,Ti} = copy(M)
SparseMatrixCSC{Tv,Ti}(M::Transpose{Tv,SparseMatrixCSC{Tv,Ti}}) where {Tv,Ti} = copy(M)
# converting from adjoint or transpose sparse matrices to sparse matrices with different eltype
SparseMatrixCSC{Tv}(M::Adjoint{<:Any,SparseMatrixCSC}) where {Tv} = SparseMatrixCSC{Tv}(copy(M))
SparseMatrixCSC{Tv}(M::Transpose{<:Any,SparseMatrixCSC}) where {Tv} = SparseMatrixCSC{Tv}(copy(M))
SparseMatrixCSC{Tv,Ti}(M::Adjoint{<:Any,SparseMatrixCSC}) where {Tv,Ti} = SparseMatrixCSC{Tv,Ti}(copy(M))
SparseMatrixCSC{Tv,Ti}(M::Transpose{<:Any,SparseMatrixCSC}) where {Tv,Ti} = SparseMatrixCSC{Tv,Ti}(copy(M))
# converting from SparseMatrixCSC to other matrix types
function Matrix(S::SparseMatrixCSC{Tv}) where Tv
# Handle cases where zero(Tv) is not defined but the array is dense.
A = length(S) == nnz(S) ? Matrix{Tv}(undef, S.m, S.n) : zeros(Tv, S.m, S.n)
for Sj in 1:S.n
for Sk in nzrange(S, Sj)
Si = S.rowval[Sk]
Sv = S.nzval[Sk]
A[Si, Sj] = Sv
end
end
return A
end
Array(S::SparseMatrixCSC) = Matrix(S)
convert(T::Type{<:SparseMatrixCSC}, m::AbstractMatrix) = m isa T ? m : T(m)
float(S::SparseMatrixCSC) = SparseMatrixCSC(S.m, S.n, copy(S.colptr), copy(S.rowval), float.(S.nzval))
complex(S::SparseMatrixCSC) = SparseMatrixCSC(S.m, S.n, copy(S.colptr), copy(S.rowval), complex(copy(S.nzval)))
"""
sparse(A)
Convert an AbstractMatrix `A` into a sparse matrix.
# Examples
```jldoctest
julia> A = Matrix(1.0I, 3, 3)
3×3 Array{Float64,2}:
1.0 0.0 0.0
0.0 1.0 0.0
0.0 0.0 1.0
julia> sparse(A)
3×3 SparseMatrixCSC{Float64,Int64} with 3 stored entries:
[1, 1] = 1.0
[2, 2] = 1.0
[3, 3] = 1.0
```
"""
sparse(A::AbstractMatrix{Tv}) where {Tv} = convert(SparseMatrixCSC{Tv,Int}, A)
sparse(S::SparseMatrixCSC) = copy(S)
sparse_IJ_sorted!(I,J,V,m,n) = sparse_IJ_sorted!(I,J,V,m,n,+)
sparse_IJ_sorted!(I,J,V::AbstractVector{Bool},m,n) = sparse_IJ_sorted!(I,J,V,m,n,|)
function sparse_IJ_sorted!(I::AbstractVector{Ti}, J::AbstractVector{Ti},
V::AbstractVector,
m::Integer, n::Integer, combine::Function) where Ti<:Integer
@assert !has_offset_axes(I, J, V)
m = m < 0 ? 0 : m
n = n < 0 ? 0 : n
if isempty(V); return spzeros(eltype(V),Ti,m,n); end
cols = zeros(Ti, n+1)
cols[1] = 1 # For cumsum purposes
cols[J[1] + 1] = 1
lastdup = 1
ndups = 0
I_lastdup = I[1]
J_lastdup = J[1]
L = length(I)
@inbounds for k=2:L
if I[k] == I_lastdup && J[k] == J_lastdup
V[lastdup] = combine(V[lastdup], V[k])
ndups += 1
else
cols[J[k] + 1] += 1
lastdup = k-ndups
I_lastdup = I[k]
J_lastdup = J[k]
if ndups != 0
I[lastdup] = I_lastdup
V[lastdup] = V[k]
end
end
end
colptr = cumsum!(similar(cols), cols)
# Allow up to 20% slack
if ndups > 0.2*L
numnz = L-ndups
deleteat!(I, (numnz+1):L)
deleteat!(V, (numnz+1):length(V))
end
return SparseMatrixCSC(m, n, colptr, I, V)
end
"""
sparse(I, J, V,[ m, n, combine])
Create a sparse matrix `S` of dimensions `m x n` such that `S[I[k], J[k]] = V[k]`. The
`combine` function is used to combine duplicates. If `m` and `n` are not specified, they
are set to `maximum(I)` and `maximum(J)` respectively. If the `combine` function is not
supplied, `combine` defaults to `+` unless the elements of `V` are Booleans in which case
`combine` defaults to `|`. All elements of `I` must satisfy `1 <= I[k] <= m`, and all
elements of `J` must satisfy `1 <= J[k] <= n`. Numerical zeros in (`I`, `J`, `V`) are
retained as structural nonzeros; to drop numerical zeros, use [`dropzeros!`](@ref).
For additional documentation and an expert driver, see `SparseArrays.sparse!`.
# Examples
```jldoctest
julia> Is = [1; 2; 3];
julia> Js = [1; 2; 3];
julia> Vs = [1; 2; 3];
julia> sparse(Is, Js, Vs)
3×3 SparseMatrixCSC{Int64,Int64} with 3 stored entries:
[1, 1] = 1
[2, 2] = 2
[3, 3] = 3
```
"""
function sparse(I::AbstractVector{Ti}, J::AbstractVector{Ti}, V::AbstractVector{Tv}, m::Integer, n::Integer, combine) where {Tv,Ti<:Integer}
@assert !has_offset_axes(I, J, V)
coolen = length(I)
if length(J) != coolen || length(V) != coolen
throw(ArgumentError(string("the first three arguments' lengths must match, ",
"length(I) (=$(length(I))) == length(J) (= $(length(J))) == length(V) (= ",
"$(length(V)))")))
end
if m == 0 || n == 0 || coolen == 0
if coolen != 0
if n == 0
throw(ArgumentError("column indices J[k] must satisfy 1 <= J[k] <= n"))
elseif m == 0
throw(ArgumentError("row indices I[k] must satisfy 1 <= I[k] <= m"))
end
end
SparseMatrixCSC(m, n, fill(one(Ti), n+1), Vector{Ti}(), Vector{Tv}())
else
# Allocate storage for CSR form
csrrowptr = Vector{Ti}(undef, m+1)
csrcolval = Vector{Ti}(undef, coolen)
csrnzval = Vector{Tv}(undef, coolen)
# Allocate storage for the CSC form's column pointers and a necessary workspace
csccolptr = Vector{Ti}(undef, n+1)
klasttouch = Vector{Ti}(undef, n)
# Allocate empty arrays for the CSC form's row and nonzero value arrays
# The parent method called below automagically resizes these arrays
cscrowval = Vector{Ti}()
cscnzval = Vector{Tv}()
sparse!(I, J, V, m, n, combine, klasttouch,
csrrowptr, csrcolval, csrnzval,
csccolptr, cscrowval, cscnzval)
end
end
sparse(I::AbstractVector, J::AbstractVector, V::AbstractVector, m::Integer, n::Integer, combine) =
sparse(AbstractVector{Int}(I), AbstractVector{Int}(J), V, m, n, combine)
"""
sparse!(I::AbstractVector{Ti}, J::AbstractVector{Ti}, V::AbstractVector{Tv},
m::Integer, n::Integer, combine, klasttouch::Vector{Ti},
csrrowptr::Vector{Ti}, csrcolval::Vector{Ti}, csrnzval::Vector{Tv},
[csccolptr::Vector{Ti}], [cscrowval::Vector{Ti}, cscnzval::Vector{Tv}] ) where {Tv,Ti<:Integer}
Parent of and expert driver for [`sparse`](@ref);
see [`sparse`](@ref) for basic usage. This method
allows the user to provide preallocated storage for `sparse`'s intermediate objects and
result as described below. This capability enables more efficient successive construction
of [`SparseMatrixCSC`](@ref)s from coordinate representations, and also enables extraction
of an unsorted-column representation of the result's transpose at no additional cost.
This method consists of three major steps: (1) Counting-sort the provided coordinate
representation into an unsorted-row CSR form including repeated entries. (2) Sweep through
the CSR form, simultaneously calculating the desired CSC form's column-pointer array,
detecting repeated entries, and repacking the CSR form with repeated entries combined;
this stage yields an unsorted-row CSR form with no repeated entries. (3) Counting-sort the
preceding CSR form into a fully-sorted CSC form with no repeated entries.
Input arrays `csrrowptr`, `csrcolval`, and `csrnzval` constitute storage for the
intermediate CSR forms and require `length(csrrowptr) >= m + 1`,
`length(csrcolval) >= length(I)`, and `length(csrnzval >= length(I))`. Input
array `klasttouch`, workspace for the second stage, requires `length(klasttouch) >= n`.
Optional input arrays `csccolptr`, `cscrowval`, and `cscnzval` constitute storage for the
returned CSC form `S`. `csccolptr` requires `length(csccolptr) >= n + 1`. If necessary,
`cscrowval` and `cscnzval` are automatically resized to satisfy
`length(cscrowval) >= nnz(S)` and `length(cscnzval) >= nnz(S)`; hence, if `nnz(S)` is
unknown at the outset, passing in empty vectors of the appropriate type (`Vector{Ti}()`
and `Vector{Tv}()` respectively) suffices, or calling the `sparse!` method
neglecting `cscrowval` and `cscnzval`.
On return, `csrrowptr`, `csrcolval`, and `csrnzval` contain an unsorted-column
representation of the result's transpose.
You may reuse the input arrays' storage (`I`, `J`, `V`) for the output arrays
(`csccolptr`, `cscrowval`, `cscnzval`). For example, you may call
`sparse!(I, J, V, csrrowptr, csrcolval, csrnzval, I, J, V)`.
For the sake of efficiency, this method performs no argument checking beyond
`1 <= I[k] <= m` and `1 <= J[k] <= n`. Use with care. Testing with `--check-bounds=yes`
is wise.
This method runs in `O(m, n, length(I))` time. The HALFPERM algorithm described in
F. Gustavson, "Two fast algorithms for sparse matrices: multiplication and permuted
transposition," ACM TOMS 4(3), 250-269 (1978) inspired this method's use of a pair of
counting sorts.
"""
function sparse!(I::AbstractVector{Ti}, J::AbstractVector{Ti},
V::AbstractVector{Tv}, m::Integer, n::Integer, combine, klasttouch::Vector{Ti},
csrrowptr::Vector{Ti}, csrcolval::Vector{Ti}, csrnzval::Vector{Tv},
csccolptr::Vector{Ti}, cscrowval::Vector{Ti}, cscnzval::Vector{Tv}) where {Tv,Ti<:Integer}
@assert !has_offset_axes(I, J, V)
# Compute the CSR form's row counts and store them shifted forward by one in csrrowptr
fill!(csrrowptr, Ti(0))
coolen = length(I)
@inbounds for k in 1:coolen
Ik = I[k]
if 1 > Ik || m < Ik
throw(ArgumentError("row indices I[k] must satisfy 1 <= I[k] <= m"))
end
csrrowptr[Ik+1] += Ti(1)
end
# Compute the CSR form's rowptrs and store them shifted forward by one in csrrowptr
countsum = Ti(1)
csrrowptr[1] = Ti(1)
@inbounds for i in 2:(m+1)
overwritten = csrrowptr[i]
csrrowptr[i] = countsum
countsum += overwritten
end
# Counting-sort the column and nonzero values from J and V into csrcolval and csrnzval
# Tracking write positions in csrrowptr corrects the row pointers
@inbounds for k in 1:coolen
Ik, Jk = I[k], J[k]
if Ti(1) > Jk || Ti(n) < Jk
throw(ArgumentError("column indices J[k] must satisfy 1 <= J[k] <= n"))
end
csrk = csrrowptr[Ik+1]
csrrowptr[Ik+1] = csrk + Ti(1)
csrcolval[csrk] = Jk
csrnzval[csrk] = V[k]
end
# This completes the unsorted-row, has-repeats CSR form's construction
# Sweep through the CSR form, simultaneously (1) calculating the CSC form's column
# counts and storing them shifted forward by one in csccolptr; (2) detecting repeated
# entries; and (3) repacking the CSR form with the repeated entries combined.
#
# Minimizing extraneous communication and nonlocality of reference, primarily by using
# only a single auxiliary array in this step, is the key to this method's performance.
fill!(csccolptr, Ti(0))
fill!(klasttouch, Ti(0))
writek = Ti(1)
newcsrrowptri = Ti(1)
origcsrrowptri = Ti(1)
origcsrrowptrip1 = csrrowptr[2]
@inbounds for i in 1:m
for readk in origcsrrowptri:(origcsrrowptrip1-Ti(1))
j = csrcolval[readk]
if klasttouch[j] < newcsrrowptri
klasttouch[j] = writek
if writek != readk
csrcolval[writek] = j
csrnzval[writek] = csrnzval[readk]
end
writek += Ti(1)
csccolptr[j+1] += Ti(1)
else
klt = klasttouch[j]
csrnzval[klt] = combine(csrnzval[klt], csrnzval[readk])
end
end
newcsrrowptri = writek
origcsrrowptri = origcsrrowptrip1
origcsrrowptrip1 != writek && (csrrowptr[i+1] = writek)
i < m && (origcsrrowptrip1 = csrrowptr[i+2])
end
# Compute the CSC form's colptrs and store them shifted forward by one in csccolptr
countsum = Ti(1)
csccolptr[1] = Ti(1)
@inbounds for j in 2:(n+1)
overwritten = csccolptr[j]
csccolptr[j] = countsum
countsum += overwritten
end
# Now knowing the CSC form's entry count, resize cscrowval and cscnzval if necessary
cscnnz = countsum - Ti(1)
length(cscrowval) < cscnnz && resize!(cscrowval, cscnnz)
length(cscnzval) < cscnnz && resize!(cscnzval, cscnnz)
# Finally counting-sort the row and nonzero values from the CSR form into cscrowval and
# cscnzval. Tracking write positions in csccolptr corrects the column pointers.
@inbounds for i in 1:m
for csrk in csrrowptr[i]:(csrrowptr[i+1]-Ti(1))
j = csrcolval[csrk]
x = csrnzval[csrk]
csck = csccolptr[j+1]
csccolptr[j+1] = csck + Ti(1)
cscrowval[csck] = i
cscnzval[csck] = x
end
end
SparseMatrixCSC(m, n, csccolptr, cscrowval, cscnzval)
end
function sparse!(I::AbstractVector{Ti}, J::AbstractVector{Ti},
V::AbstractVector{Tv}, m::Integer, n::Integer, combine, klasttouch::Vector{Ti},
csrrowptr::Vector{Ti}, csrcolval::Vector{Ti}, csrnzval::Vector{Tv},
csccolptr::Vector{Ti}) where {Tv,Ti<:Integer}
sparse!(I, J, V, m, n, combine, klasttouch,
csrrowptr, csrcolval, csrnzval,
csccolptr, Vector{Ti}(), Vector{Tv}())
end
function sparse!(I::AbstractVector{Ti}, J::AbstractVector{Ti},
V::AbstractVector{Tv}, m::Integer, n::Integer, combine, klasttouch::Vector{Ti},
csrrowptr::Vector{Ti}, csrcolval::Vector{Ti}, csrnzval::Vector{Tv}) where {Tv,Ti<:Integer}
sparse!(I, J, V, m, n, combine, klasttouch,
csrrowptr, csrcolval, csrnzval,
Vector{Ti}(undef, n+1), Vector{Ti}(), Vector{Tv}())
end
dimlub(I) = isempty(I) ? 0 : Int(maximum(I)) #least upper bound on required sparse matrix dimension
sparse(I,J,v::Number) = sparse(I, J, fill(v,length(I)))
sparse(I,J,V::AbstractVector) = sparse(I, J, V, dimlub(I), dimlub(J))
sparse(I,J,v::Number,m,n) = sparse(I, J, fill(v,length(I)), Int(m), Int(n))
sparse(I,J,V::AbstractVector,m,n) = sparse(I, J, V, Int(m), Int(n), +)
sparse(I,J,V::AbstractVector{Bool},m,n) = sparse(I, J, V, Int(m), Int(n), |)
sparse(I,J,v::Number,m,n,combine::Function) = sparse(I, J, fill(v,length(I)), Int(m), Int(n), combine)
function sparse(T::SymTridiagonal)
m = length(T.dv)
return sparse([1:m;2:m;1:m-1],[1:m;1:m-1;2:m],[T.dv;T.ev;T.ev], Int(m), Int(m))
end
function sparse(T::Tridiagonal)
m = length(T.d)
return sparse([1:m;2:m;1:m-1],[1:m;1:m-1;2:m],[T.d;T.dl;T.du], Int(m), Int(m))
end
function sparse(B::Bidiagonal)
m = length(B.dv)
B.uplo == 'U' || return sparse([1:m;2:m],[1:m;1:m-1],[B.dv;B.ev], Int(m), Int(m)) # lower bidiagonal
return sparse([1:m;1:m-1],[1:m;2:m],[B.dv;B.ev], Int(m), Int(m)) # upper bidiagonal
end
function sparse(D::Diagonal{T}) where T
m = length(D.diag)
return SparseMatrixCSC(m, m, Vector(1:(m+1)), Vector(1:m), Vector{T}(D.diag))
end
## Transposition and permutation methods
"""
halfperm!(X::SparseMatrixCSC{Tv,Ti}, A::SparseMatrixCSC{Tv,Ti},
q::AbstractVector{<:Integer}, f::Function = identity) where {Tv,Ti}
Column-permute and transpose `A`, simultaneously applying `f` to each entry of `A`, storing
the result `(f(A)Q)^T` (`map(f, transpose(A[:,q]))`) in `X`.
`X`'s dimensions must match those of `transpose(A)` (`X.m == A.n` and `X.n == A.m`), and `X`
must have enough storage to accommodate all allocated entries in `A` (`length(X.rowval) >= nnz(A)`
and `length(X.nzval) >= nnz(A)`). Column-permutation `q`'s length must match `A`'s column
count (`length(q) == A.n`).
This method is the parent of several methods performing transposition and permutation
operations on [`SparseMatrixCSC`](@ref)s. As this method performs no argument checking,
prefer the safer child methods (`[c]transpose[!]`, `permute[!]`) to direct use.
This method implements the `HALFPERM` algorithm described in F. Gustavson, "Two fast
algorithms for sparse matrices: multiplication and permuted transposition," ACM TOMS 4(3),
250-269 (1978). The algorithm runs in `O(A.m, A.n, nnz(A))` time and requires no space
beyond that passed in.
"""
function halfperm!(X::SparseMatrixCSC{Tv,Ti}, A::SparseMatrixCSC{Tv,Ti},
q::AbstractVector{<:Integer}, f::Function = identity) where {Tv,Ti}
_computecolptrs_halfperm!(X, A)
_distributevals_halfperm!(X, A, q, f)
return X
end
"""
Helper method for `halfperm!`. Computes `transpose(A[:,q])`'s column pointers, storing them
shifted one position forward in `X.colptr`; `_distributevals_halfperm!` fixes this shift.
"""
function _computecolptrs_halfperm!(X::SparseMatrixCSC{Tv,Ti}, A::SparseMatrixCSC{Tv,Ti}) where {Tv,Ti}
# Compute `transpose(A[:,q])`'s column counts. Store shifted forward one position in X.colptr.
fill!(X.colptr, 0)
@inbounds for k in 1:nnz(A)
X.colptr[A.rowval[k] + 1] += 1
end
# Compute `transpose(A[:,q])`'s column pointers. Store shifted forward one position in X.colptr.
X.colptr[1] = 1
countsum = 1
@inbounds for k in 2:(A.m + 1)
overwritten = X.colptr[k]
X.colptr[k] = countsum
countsum += overwritten
end
end
"""
Helper method for `halfperm!`. With `transpose(A[:,q])`'s column pointers shifted one
position forward in `X.colptr`, computes `map(f, transpose(A[:,q]))` by appropriately
distributing `A.rowval` and `f`-transformed `A.nzval` into `X.rowval` and `X.nzval`
respectively. Simultaneously fixes the one-position-forward shift in `X.colptr`.
"""
@noinline function _distributevals_halfperm!(X::SparseMatrixCSC{Tv,Ti},
A::SparseMatrixCSC{Tv,Ti}, q::AbstractVector{<:Integer}, f::Function) where {Tv,Ti}
@inbounds for Xi in 1:A.n
Aj = q[Xi]
for Ak in nzrange(A, Aj)
Ai = A.rowval[Ak]
Xk = X.colptr[Ai + 1]
X.rowval[Xk] = Xi
X.nzval[Xk] = f(A.nzval[Ak])
X.colptr[Ai + 1] += 1
end
end
return # kill potential type instability
end
function ftranspose!(X::SparseMatrixCSC{Tv,Ti}, A::SparseMatrixCSC{Tv,Ti}, f::Function) where {Tv,Ti}
# Check compatibility of source argument A and destination argument X
if X.n != A.m
throw(DimensionMismatch(string("destination argument `X`'s column count, ",
"`X.n (= $(X.n))`, must match source argument `A`'s row count, `A.m (= $(A.m))`")))
elseif X.m != A.n
throw(DimensionMismatch(string("destination argument `X`'s row count,
`X.m (= $(X.m))`, must match source argument `A`'s column count, `A.n (= $(A.n))`")))
elseif length(X.rowval) < nnz(A)
throw(ArgumentError(string("the length of destination argument `X`'s `rowval` ",
"array, `length(X.rowval) (= $(length(X.rowval)))`, must be greater than or ",
"equal to source argument `A`'s allocated entry count, `nnz(A) (= $(nnz(A)))`")))
elseif length(X.nzval) < nnz(A)
throw(ArgumentError(string("the length of destination argument `X`'s `nzval` ",
"array, `length(X.nzval) (= $(length(X.nzval)))`, must be greater than or ",
"equal to source argument `A`'s allocated entry count, `nnz(A) (= $(nnz(A)))`")))
end
halfperm!(X, A, 1:A.n, f)
end
transpose!(X::SparseMatrixCSC{Tv,Ti}, A::SparseMatrixCSC{Tv,Ti}) where {Tv,Ti} = ftranspose!(X, A, identity)
adjoint!(X::SparseMatrixCSC{Tv,Ti}, A::SparseMatrixCSC{Tv,Ti}) where {Tv,Ti} = ftranspose!(X, A, conj)
function ftranspose(A::SparseMatrixCSC{Tv,Ti}, f::Function) where {Tv,Ti}
X = SparseMatrixCSC(A.n, A.m,
Vector{Ti}(undef, A.m+1),
Vector{Ti}(undef, nnz(A)),
Vector{Tv}(undef, nnz(A)))
halfperm!(X, A, 1:A.n, f)
end
adjoint(A::SparseMatrixCSC) = Adjoint(A)
transpose(A::SparseMatrixCSC) = Transpose(A)
Base.copy(A::Adjoint{<:Any,<:SparseMatrixCSC}) = ftranspose(A.parent, x -> copy(adjoint(x)))
Base.copy(A::Transpose{<:Any,<:SparseMatrixCSC}) = ftranspose(A.parent, x -> copy(transpose(x)))
function Base.permutedims(A::SparseMatrixCSC, (a,b))
(a, b) == (2, 1) && return ftranspose(A, identity)
(a, b) == (1, 2) && return copy(A)
throw(ArgumentError("no valid permutation of dimensions"))
end
"""
unchecked_noalias_permute!(X::SparseMatrixCSC{Tv,Ti},
A::SparseMatrixCSC{Tv,Ti}, p::AbstractVector{<:Integer},
q::AbstractVector{<:Integer}, C::SparseMatrixCSC{Tv,Ti}) where {Tv,Ti}
See [`permute!`](@ref) for basic usage. Parent of `permute[!]`
methods operating on `SparseMatrixCSC`s that assume none of `X`, `A`, and `C` alias each
other. As this method performs no argument checking, prefer the safer child methods
(`permute[!]`) to direct use.
This method consists of two major steps: (1) Column-permute (`Q`,`I[:,q]`) and transpose `A`
to generate intermediate result `(AQ)^T` (`transpose(A[:,q])`) in `C`. (2) Column-permute
(`P^T`, I[:,p]) and transpose intermediate result `(AQ)^T` to generate result
`((AQ)^T P^T)^T = PAQ` (`A[p,q]`) in `X`.
The first step is a call to `halfperm!`, and the second is a variant on `halfperm!` that
avoids an unnecessary length-`nnz(A)` array-sweep and associated recomputation of column
pointers. See [`halfperm!`](:func:SparseArrays.halfperm!) for additional algorithmic
information.
See also: `unchecked_aliasing_permute!`
"""
function unchecked_noalias_permute!(X::SparseMatrixCSC{Tv,Ti},
A::SparseMatrixCSC{Tv,Ti}, p::AbstractVector{<:Integer},
q::AbstractVector{<:Integer}, C::SparseMatrixCSC{Tv,Ti}) where {Tv,Ti}
halfperm!(C, A, q)
_computecolptrs_permute!(X, A, q, X.colptr)
_distributevals_halfperm!(X, C, p, identity)
return X
end
"""
unchecked_aliasing_permute!(A::SparseMatrixCSC{Tv,Ti},
p::AbstractVector{<:Integer}, q::AbstractVector{<:Integer},
C::SparseMatrixCSC{Tv,Ti}, workcolptr::Vector{Ti}) where {Tv,Ti}
See [`permute!`](@ref) for basic usage. Parent of `permute!`
methods operating on [`SparseMatrixCSC`](@ref)s where the source and destination matrices
are the same. See `unchecked_noalias_permute!`
for additional information; these methods are identical but for this method's requirement of
the additional `workcolptr`, `length(workcolptr) >= A.n + 1`, which enables efficient
handling of the source-destination aliasing.
"""
function unchecked_aliasing_permute!(A::SparseMatrixCSC{Tv,Ti},
p::AbstractVector{<:Integer}, q::AbstractVector{<:Integer},
C::SparseMatrixCSC{Tv,Ti}, workcolptr::Vector{Ti}) where {Tv,Ti}
halfperm!(C, A, q)
_computecolptrs_permute!(A, A, q, workcolptr)
_distributevals_halfperm!(A, C, p, identity)
return A
end
"""
Helper method for `unchecked_noalias_permute!` and `unchecked_aliasing_permute!`.
Computes `PAQ`'s column pointers, storing them shifted one position forward in `X.colptr`;
`_distributevals_halfperm!` fixes this shift. Saves some work relative to
`_computecolptrs_halfperm!` as described in `uncheckednoalias_permute!`'s documentation.
"""
function _computecolptrs_permute!(X::SparseMatrixCSC{Tv,Ti},
A::SparseMatrixCSC{Tv,Ti}, q::AbstractVector{<:Integer}, workcolptr::Vector{Ti}) where {Tv,Ti}
# Compute `A[p,q]`'s column counts. Store shifted forward one position in workcolptr.
@inbounds for k in 1:A.n
workcolptr[k+1] = A.colptr[q[k] + 1] - A.colptr[q[k]]
end
# Compute `A[p,q]`'s column pointers. Store shifted forward one position in X.colptr.
X.colptr[1] = 1
countsum = 1
@inbounds for k in 2:(X.n + 1)
overwritten = workcolptr[k]
X.colptr[k] = countsum
countsum += overwritten
end
end
"""
Helper method for `permute` and `permute!` methods operating on `SparseMatrixCSC`s.
Checks compatibility of source argument `A`, row-permutation argument `p`, and
column-permutation argument `q`.
"""
function _checkargs_sourcecompatperms_permute!(A::SparseMatrixCSC,
p::AbstractVector{<:Integer}, q::AbstractVector{<:Integer})
@assert !has_offset_axes(p, q)
if length(q) != A.n
throw(DimensionMismatch(string("the length of column-permutation argument `q`, ",
"`length(q) (= $(length(q)))`, must match source argument `A`'s column ",
"count, `A.n (= $(A.n))`")))
elseif length(p) != A.m
throw(DimensionMismatch(string("the length of row-permutation argument `p`, ",
"`length(p) (= $(length(p)))`, must match source argument `A`'s row count, ",
"`A.m (= $(A.m))`")))
end
end
"""
Helper method for `permute` and `permute!` methods operating on `SparseMatrixCSC`s.
Checks whether row- and column- permutation arguments `p` and `q` are valid permutations.
"""
function _checkargs_permutationsvalid_permute!(
p::AbstractVector{<:Integer}, pcheckspace::Vector{Ti},
q::AbstractVector{<:Integer}, qcheckspace::Vector{Ti}) where Ti<:Integer
if !_ispermutationvalid_permute!(p, pcheckspace)
throw(ArgumentError("row-permutation argument `p` must be a valid permutation"))
elseif !_ispermutationvalid_permute!(q, qcheckspace)
throw(ArgumentError("column-permutation argument `q` must be a valid permutation"))
end
end
function _ispermutationvalid_permute!(perm::AbstractVector{<:Integer},
checkspace::Vector{<:Integer})
@assert !has_offset_axes(perm)
n = length(perm)
checkspace[1:n] .= 0
for k in perm
(0 < k ≤ n) && ((checkspace[k] ⊻= 1) == 1) || return false
end
return true