Skip to content
This repository was archived by the owner on Aug 22, 2025. It is now read-only.

Commit 9a2d26f

Browse files
committed
Merge branch 'master' of https://github.com/JuliaDiffEq/SparseDiffTools.jl into coloring
2 parents 495916b + 05a90f9 commit 9a2d26f

File tree

2 files changed

+4
-2
lines changed

2 files changed

+4
-2
lines changed

src/coloring/high_level.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -74,9 +74,9 @@ function matrix_colors(A::BandedBlockBandedMatrix)
7474
cols=[blocksize(A,(1,i))[2] for i in 1:nblock]
7575
blockcolors=_cycle(1:blockwidth,nblock)
7676
#the reserved number of colors of a block is the min of subblockwidth and the largest length of columns of blocks with the same block color
77-
ncolors=[min(subblockwidth,maximum(cols[i:blockwidth:nblock])) for i in 1:blockwidth]
77+
ncolors=[min(subblockwidth,maximum(cols[i:blockwidth:nblock])) for i in 1:min(blockwidth,nblock)]
7878
endinds=cumsum(ncolors)
79-
startinds=[endinds[i]-ncolors[i]+1 for i in 1:blockwidth]
79+
startinds=[endinds[i]-ncolors[i]+1 for i in 1:min(blockwidth,nblock)]
8080
colors=[_cycle(startinds[blockcolors[i]]:endinds[blockcolors[i]],cols[i]) for i in 1:nblock]
8181
vcat(colors...)
8282
end

test/test_specialmatrices.jl

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -19,6 +19,7 @@ blockbanded1=BlockBandedMatrix(dense,([1,2,3,4],[4,3,2,1]),(1,0))
1919
blockbanded2=BlockBandedMatrix(dense,([4,3,2,1],[1,2,3,4]),(1,1))
2020
bandedblockbanded1=BandedBlockBandedMatrix(dense,([1,2,3,4],[4,3,2,1]),(1,0),(1,1))
2121
bandedblockbanded2=BandedBlockBandedMatrix(dense,([4,3,2,1],[1,2,3,4]),(1,1),(1,0))
22+
bandedblockbanded3 = BandedBlockBandedMatrix(Zeros(2 * 4 , 2 * 4), ([4, 4] ,[4, 4]), (1, 1), (1, 1))
2223

2324
@test matrix_colors(dense)==1:n
2425
@test matrix_colors(uptri)==1:n
@@ -35,6 +36,7 @@ bandedblockbanded2=BandedBlockBandedMatrix(dense,([4,3,2,1],[1,2,3,4]),(1,1),(1,
3536
@test matrix_colors(blockbanded2)==[1,5,6,7,8,9,1,2,3,4]
3637
@test matrix_colors(bandedblockbanded1)==[1,2,3,1,4,5,6,1,2,4]
3738
@test matrix_colors(bandedblockbanded2)==[1,3,4,5,6,5,1,2,1,2]
39+
@test matrix_colors(bandedblockbanded3)==[1,2,3,1,4,5,6,4]
3840

3941
function _testvalidity(A)
4042
colorvec=matrix_colors(A)

0 commit comments

Comments
 (0)