Internals
The docstrings on this page describe internals, they are not part of the public API.
Graph storage
SparseMatrixColorings.SparsityPatternCSC
— TypeSparsityPatternCSC{Ti<:Integer}
Store a sparse matrix (in CSC) without its values, keeping only the pattern of nonzeros.
Fields
Copied from SparseMatrixCSC
:
m::Int
: number of rowsn::Int
: number of columnscolptr::Vector{Ti}
: columnj
is incolptr[j]:(colptr[j+1]-1)
rowval::Vector{Ti}
: row indices of stored values
SparseMatrixColorings.AdjacencyGraph
— TypeAdjacencyGraph{T}
Undirected graph without self-loops representing the nonzeros of a symmetric matrix (typically a Hessian matrix).
The adjacency graph of a symmetrix matric A ∈ ℝ^{n × n}
is G(A) = (V, E)
where
V = 1:n
is the set of rows or columnsi
/j
(i, j) ∈ E
wheneverA[i, j] ≠ 0
andi ≠ j
Constructors
AdjacencyGraph(A::SparseMatrixCSC)
Fields
S::SparsityPatternCSC{T}
References
What Color Is Your Jacobian? SparsityPatternCSC Coloring for Computing Derivatives, Gebremedhin et al. (2005)
SparseMatrixColorings.BipartiteGraph
— TypeBipartiteGraph{T}
Undirected bipartite graph representing the nonzeros of a non-symmetric matrix (typically a Jacobian matrix).
The bipartite graph of a matrix A ∈ ℝ^{m × n}
is Gb(A) = (V₁, V₂, E)
where
V₁ = 1:m
is the set of rowsi
V₂ = 1:n
is the set of columnsj
(i, j) ∈ E
wheneverA[i, j] ≠ 0
A BipartiteGraph
has two sets of vertices, one for the rows of A
(which we call side 1
) and one for the columns (which we call side 2
).
Constructors
BipartiteGraph(A::SparseMatrixCSC; symmetric_pattern=false)
When symmetric_pattern
is true
, this construction is more efficient.
Fields
S1::SparsityPatternCSC{T}
: maps vertices on side1
to their neighborsS2::SparsityPatternCSC{T}
: maps vertices on side2
to their neighbors
References
What Color Is Your Jacobian? SparsityPatternCSC Coloring for Computing Derivatives, Gebremedhin et al. (2005)
SparseMatrixColorings.vertices
— Functionvertices(bg::BipartiteGraph, Val(side))
Return the list of vertices of bg
from the specified side
as a range 1:n
.
SparseMatrixColorings.neighbors
— Functionneighbors(bg::BipartiteGraph, Val(side), v::Integer)
Return the neighbors of v
(a vertex from the specified side
, 1
or 2
), in the graph bg
.
Base.transpose
— Functiontranspose(S::SparsityPatternCSC)
Return a SparsityPatternCSC
corresponding to the transpose of S
.
Low-level coloring
SparseMatrixColorings.partial_distance2_coloring
— Functionpartial_distance2_coloring(bg::BipartiteGraph, ::Val{side}, order::AbstractOrder)
Compute a distance-2 coloring of the given side
(1
or 2
) in the bipartite graph bg
and return a vector of integer colors.
A distance-2 coloring is such that two vertices have different colors if they are at distance at most 2.
The vertices are colored in a greedy fashion, following the order
supplied.
See also
References
What Color Is Your Jacobian? Graph Coloring for Computing Derivatives, Gebremedhin et al. (2005), Algorithm 3.2
SparseMatrixColorings.symmetric_coefficient
— Functionsymmetric_coefficient(
i::Integer, j::Integer,
color::AbstractVector{<:Integer},
star_set::StarSet
)
Return the indices (k, c)
such that A[i, j] = B[k, c]
, where A
is the uncompressed symmetric matrix and B
is the column-compressed matrix.
This function corresponds to algorithm DirectRecover2
in the paper.
References
Efficient Computation of Sparse Hessians Using Coloring and Automatic Differentiation, Gebremedhin et al. (2009), Figure 3
SparseMatrixColorings.star_coloring
— Functionstar_coloring(g::AdjacencyGraph, order::AbstractOrder)
Compute a star coloring of all vertices in the adjacency graph g
and return a tuple (color, star_set)
, where
color
is the vector of integer colorsstar_set
is aStarSet
encoding the set of 2-colored stars
A star coloring is a distance-1 coloring such that every path on 4 vertices uses at least 3 colors.
The vertices are colored in a greedy fashion, following the order
supplied.
See also
References
New Acyclic and Star Coloring Algorithms with Application to Computing Hessians, Gebremedhin et al. (2007), Algorithm 4.1
SparseMatrixColorings.acyclic_coloring
— Functionacyclic_coloring(g::AdjacencyGraph, order::AbstractOrder)
Compute an acyclic coloring of all vertices in the adjacency graph g
and return a tuple (color, tree_set)
, where
color
is the vector of integer colorstree_set
is aTreeSet
encoding the set of 2-colored trees
An acyclic coloring is a distance-1 coloring with the further restriction that every cycle uses at least 3 colors.
The vertices are colored in a greedy fashion, following the order
supplied.
See also
References
New Acyclic and Star Coloring Algorithms with Application to Computing Hessians, Gebremedhin et al. (2007), Algorithm 3.1
SparseMatrixColorings.group_by_color
— Functiongroup_by_color(color::Vector{Int})
Create a color-indexed vector group
such that i ∈ group[c]
iff color[i] == c
.
Assumes the colors are contiguously numbered from 1
to some cmax
.
SparseMatrixColorings.StarSet
— TypeStarSet
Encode a set of 2-colored stars resulting from the star_coloring
algorithm.
Fields
star::Dict{Tuple{Int64, Int64}, Int64}
: a mapping from edges (pair of vertices) to their star indexhub::Vector{Int64}
: a mapping from star indices to their hub (undefined hubs for single-edge stars are the negative value of one of the vertices, picked arbitrarily)spokes::Vector{Vector{Int64}}
: a mapping from star indices to the vector of their spokes
SparseMatrixColorings.TreeSet
— TypeTreeSet
Encode a set of 2-colored trees resulting from the acyclic_coloring
algorithm.
Fields
forest::DataStructures.DisjointSets{Tuple{Int64, Int64}}
: a forest of two-colored trees
Concrete coloring results
SparseMatrixColorings.ColumnColoringResult
— Typestruct ColumnColoringResult{M<:(AbstractMatrix), G<:SparseMatrixColorings.BipartiteGraph, V} <: AbstractColoringResult{:nonsymmetric, :column, :direct}
Storage for the result of a column coloring with direct decompression.
Fields
A::AbstractMatrix
: matrix that was coloredbg::SparseMatrixColorings.BipartiteGraph
: bipartite graph that was used for coloringcolor::Vector{Int64}
: one integer color for each column or row (depending onpartition
)group::Any
: color groups for columns or rows (depending onpartition
)compressed_indices::Vector{Int64}
: flattened indices mapping the compressed matrixB
to the uncompressed matrixA
whenA isa SparseMatrixCSC
. They satisfynonzeros(A)[k] = vec(B)[compressed_indices[k]]
See also
SparseMatrixColorings.RowColoringResult
— Typestruct RowColoringResult{M<:(AbstractMatrix), G<:SparseMatrixColorings.BipartiteGraph, V} <: AbstractColoringResult{:nonsymmetric, :row, :direct}
Storage for the result of a row coloring with direct decompression.
Fields
See the docstring of ColumnColoringResult
.
A::AbstractMatrix
bg::SparseMatrixColorings.BipartiteGraph
color::Vector{Int64}
group::Any
compressed_indices::Vector{Int64}
See also
SparseMatrixColorings.StarSetColoringResult
— Typestruct StarSetColoringResult{M<:(AbstractMatrix), G<:SparseMatrixColorings.AdjacencyGraph, V} <: AbstractColoringResult{:symmetric, :column, :direct}
Storage for the result of a symmetric coloring with direct decompression.
Fields
See the docstring of ColumnColoringResult
.
A::AbstractMatrix
ag::SparseMatrixColorings.AdjacencyGraph
color::Vector{Int64}
group::Any
star_set::SparseMatrixColorings.StarSet
compressed_indices::Vector{Int64}
See also
SparseMatrixColorings.TreeSetColoringResult
— Typestruct TreeSetColoringResult{M<:(AbstractMatrix), G<:SparseMatrixColorings.AdjacencyGraph, V, R} <: AbstractColoringResult{:symmetric, :column, :substitution}
Storage for the result of a symmetric coloring with decompression by substitution.
Fields
See the docstring of ColumnColoringResult
.
A::AbstractMatrix
ag::SparseMatrixColorings.AdjacencyGraph
color::Vector{Int64}
group::Any
vertices_by_tree::Vector{Vector{Int64}}
reverse_bfs_orders::Vector{Vector{Tuple{Int64, Int64}}}
diagonal_indices::Vector{Int64}
diagonal_nzind::Vector{Int64}
lower_triangle_offsets::Vector{Int64}
upper_triangle_offsets::Vector{Int64}
buffer::Vector
See also
SparseMatrixColorings.LinearSystemColoringResult
— Typestruct LinearSystemColoringResult{M<:(AbstractMatrix), G<:SparseMatrixColorings.AdjacencyGraph, V, R, F} <: AbstractColoringResult{:symmetric, :column, :substitution}
Storage for the result of a symmetric coloring with any decompression.
Fields
See the docstring of ColumnColoringResult
.
A::AbstractMatrix
ag::SparseMatrixColorings.AdjacencyGraph
color::Vector{Int64}
group::Any
strict_upper_nonzero_inds::Vector{Tuple{Int64, Int64}}
strict_upper_nonzeros_A::Vector
T_factorization::Any
See also
SparseMatrixColorings.BicoloringResult
— Typestruct BicoloringResult{M<:(AbstractMatrix), G<:SparseMatrixColorings.AdjacencyGraph, decompression, V, SR<:AbstractColoringResult{:symmetric, :column, decompression}, R} <: AbstractColoringResult{:nonsymmetric, :bidirectional, decompression}
Storage for the result of a bidirectional coloring with direct or substitution decompression, based on the symmetric coloring of a 2x2 block matrix.
Fields
A::AbstractMatrix
: matrix that was coloredabg::SparseMatrixColorings.AdjacencyGraph
: adjacency graph that was used for coloring (constructed from the bipartite graph)column_color::Vector{Int64}
: one integer color for each columnrow_color::Vector{Int64}
: one integer color for each rowcolumn_group::Any
: color groups for columnsrow_group::Any
: color groups for rowssymmetric_result::AbstractColoringResult{:symmetric, :column}
: result for the coloring of the symmetric 2x2 block matrixcol_color_ind::Dict{Int64, Int64}
: column color to indexrow_color_ind::Dict{Int64, Int64}
: row color to indexBr_and_Bc::Matrix
: combination ofBr
andBc
(almost a concatenation up to color remapping)large_colptr::Vector{Int64}
: CSC storage ofA_and_noAᵀ -
colptr`large_rowval::Vector{Int64}
: CSC storage ofA_and_noAᵀ -
rowval`
See also
SparseMatrixColorings.remap_colors
— Functionremap_colors(color::Vector{Int})
Renumber the colors in color
using their index in the vector sort(unique(color))
, so that they are forced to go from 1
to some cmax
contiguously.
Return a tuple (remapped_colors, color_to_ind)
such that remapped_colors
is a vector containing the renumbered colors and color_to_ind
is a dictionary giving the translation between old and new color numberings.
For all vertex indices i
we have:
remapped_color[i] = color_to_ind[color[i]]
Testing
SparseMatrixColorings.directly_recoverable_columns
— Functiondirectly_recoverable_columns(
A::AbstractMatrix, color::AbstractVector{<:Integer}
verbose=false
)
Return true
if coloring the columns of the symmetric matrix A
with the vector color
results in a column-compressed representation that preserves every unique value, thus making direct recovery possible.
This function is not coded with efficiency in mind, it is designed for small-scale tests.
References
What Color Is Your Jacobian? Graph Coloring for Computing Derivatives, Gebremedhin et al. (2005)
SparseMatrixColorings.symmetrically_orthogonal_columns
— Functionsymmetrically_orthogonal_columns(
A::AbstractMatrix, color::AbstractVector{<:Integer};
verbose=false
)
Return true
if coloring the columns of the symmetric matrix A
with the vector color
results in a partition that is symmetrically orthogonal, and false
otherwise.
A partition of the columns of a symmetrix matrix A
is symmetrically orthogonal if, for every nonzero element A[i, j]
, either of the following statements holds:
- the group containing the column
A[:, j]
has no other column with a nonzero in rowi
- the group containing the column
A[:, i]
has no other column with a nonzero in rowj
This function is not coded with efficiency in mind, it is designed for small-scale tests.
References
What Color Is Your Jacobian? Graph Coloring for Computing Derivatives, Gebremedhin et al. (2005)
SparseMatrixColorings.structurally_orthogonal_columns
— Functionstructurally_orthogonal_columns(
A::AbstractMatrix, color::AbstractVector{<:Integer}
verbose=false
)
Return true
if coloring the columns of the matrix A
with the vector color
results in a partition that is structurally orthogonal, and false
otherwise.
A partition of the columns of a matrix A
is structurally orthogonal if, for every nonzero element A[i, j]
, the group containing column A[:, j]
has no other column with a nonzero in row i
.
This function is not coded with efficiency in mind, it is designed for small-scale tests.
References
What Color Is Your Jacobian? Graph Coloring for Computing Derivatives, Gebremedhin et al. (2005)
SparseMatrixColorings.valid_dynamic_order
— Functionvalid_dynamic_order(g::AdjacencyGraph, π::AbstractVector{Int}, order::DynamicDegreeBasedOrder)
valid_dynamic_order(bg::AdjacencyGraph, ::Val{side}, π::AbstractVector{Int}, order::DynamicDegreeBasedOrder)
Check that a permutation π
corresponds to a valid application of a DynamicDegreeBasedOrder
.
This is done by checking, for each ordered vertex, that its back- or forward-degree was the smallest or largest among the remaining vertices (the specifics depend on the order parameters).
This function is not coded with efficiency in mind, it is designed for small-scale tests.
Matrix handling
SparseMatrixColorings.respectful_similar
— Functionrespectful_similar(A::AbstractMatrix)
respectful_similar(A::AbstractMatrix, ::Type{T})
Like Base.similar
but returns a transpose or adjoint when A
is a transpose or adjoint.
SparseMatrixColorings.matrix_versions
— Functionmatrix_versions(A::AbstractMatrix)
Return various versions of the same matrix:
- dense and sparse
- transpose and adjoint
Used for internal testing.
SparseMatrixColorings.same_pattern
— Functionsame_pattern(A, B)
Perform a partial equality check on the sparsity patterns of A
and B
:
- if the return is
true
, they might have the same sparsity pattern but we're not sure - if the return is
false
, they definitely don't have the same sparsity pattern
Visualization
SparseMatrixColorings.show_colors
— Functionshow_colors(result; kwargs...)
Create a visualization for an AbstractColoringResult
, with the help of the the JuliaImages ecosystem.
- For
:column
or:row
colorings, it returns a tuple(A_img, B_img)
. - For
:bidirectional
colorings, it returns a tuple(A_img, Br_img, Bc_img)
.
This function is implemented in a package extension, using it requires loading Colors.jl.
Keyword arguments
colorscheme
: colors used for non-zero matrix entries. This can be a vector ofColorant
s or a subsampled scheme from ColorSchemes.jl.background::Colorant
: color used for zero matrix entries and pad. Defaults toRGBA(0,0,0,0)
, a transparent background.scale::Int
: scale the size of matrix entries toscale × scale
pixels. Defaults to1
.pad::Int
: set padding between matrix entries, in pixels. Defaults to0
.
For a matrix of size (n, m)
, the resulting output will be of size (n * (scale + pad) + pad, m * (scale + pad) + pad)
.
Examples
SparseMatrixColorings.Example
— Typestruct Example{TA<:(AbstractMatrix), TB<:(AbstractMatrix)}
Example coloring problem from one of our reference articles.
Used for internal testing.
Fields
A::AbstractMatrix
: decompressed matrixB::AbstractMatrix
: column-compressed matrixcolor::Vector{Int64}
: vector of colors
SparseMatrixColorings.what_fig_41
— Functionwhat_fig_41()
Construct an Example
from Figure 4.1 of "What color is your Jacobian?", where the nonzero entries are filled with unique values.
SparseMatrixColorings.what_fig_61
— Functionwhat_fig_61()
Construct an Example
from Figure 6.1 of "What color is your Jacobian?", where the nonzero entries are filled with unique values.
SparseMatrixColorings.efficient_fig_1
— Functionefficient_fig_1()
Construct an Example
from Figure 1 of "Efficient computation of sparse hessians using coloring and AD", where the nonzero entries are filled with unique values.
SparseMatrixColorings.efficient_fig_4
— Functionefficient_fig_4()
Construct an Example
from Figure 4 of "Efficient computation of sparse hessians using coloring and AD", where the nonzero entries are filled with unique values.