API
SparseDiffTools.ApproximateJacobianSparsity
— TypeApproximateJacobianSparsity(; ntrials = 5, rng = Random.default_rng(),
epsilon = nothing, alg = GreedyD1Color())
Use ntrials
random vectors to compute the sparsity pattern of the Jacobian. This is an approximate method and the sparsity pattern may not be exact.
Keyword Arguments
- `ntrials`: The number of random vectors to use for computing the sparsity pattern
- `rng`: The random number generator used for generating the random vectors
- `alg`: The algorithm used for computing the matrix colors
- `epsilon`: For Finite Differencing based Jacobian Approximations, any number smaller
than `epsilon` is considered to be zero. If `nothing` is specified, then this value
is calculated as `100 * eps(eltype(x))`
SparseDiffTools.JacPrototypeSparsityDetection
— TypeJacPrototypeSparsityDetection(; jac_prototype, alg = GreedyD1Color())
Use a pre-specified jac_prototype
to compute the matrix colors of the Jacobian.
Keyword Arguments
- `jac_prototype`: The prototype Jacobian used for computing the matrix colors
- `alg`: The algorithm used for computing the matrix colors
See Also: SymbolicsSparsityDetection, PrecomputedJacobianColorvec
SparseDiffTools.PrecomputedJacobianColorvec
— TypePrecomputedJacobianColorvec(jac_prototype, row_colorvec, col_colorvec)
Use a pre-specified colorvec
which can be directly used for sparse differentiation. Based on whether a reverse mode or forward mode or finite differences is used, the corresponding row_colorvec
or col_colorvec
is used. Atmost one of them can be set to nothing
.
Arguments
- `jac_prototype`: The prototype Jacobian used for computing structural nonzeros
- `row_colorvec`: The row colorvec of the Jacobian
- `col_colorvec`: The column colorvec of the Jacobian
See Also: SymbolicsSparsityDetection, JacPrototypeSparsityDetection
SparseDiffTools.PrecomputedJacobianColorvec
— MethodPrecomputedJacobianColorvec(; jac_prototype, partition_by_rows::Bool = false,
colorvec = missing, row_colorvec = missing, col_colorvec = missing)
Use a pre-specified colorvec
which can be directly used for sparse differentiation. Based on whether a reverse mode or forward mode or finite differences is used, the corresponding row_colorvec
or col_colorvec
is used. Atmost one of them can be set to nothing
.
Keyword Arguments
- `jac_prototype`: The prototype Jacobian used for computing structural nonzeros
- `partition_by_rows`: Whether to partition the Jacobian by rows or columns (row
partitioning is used for reverse mode AD)
- `colorvec`: The colorvec of the Jacobian. If `partition_by_rows` is `true` then this
is the row colorvec, otherwise it is the column colorvec
- `row_colorvec`: The row colorvec of the Jacobian
- `col_colorvec`: The column colorvec of the Jacobian
See Also: SymbolicsSparsityDetection, JacPrototypeSparsityDetection
SparseDiffTools.SymbolicsSparsityDetection
— TypeSymbolicsSparsityDetection(; alg = GreedyD1Color())
Use Symbolics to compute the sparsity pattern of the Jacobian. This requires Symbolics.jl
to be explicitly loaded.
Keyword Arguments
- `alg`: The algorithm used for computing the matrix colors
See Also: JacPrototypeSparsityDetection, PrecomputedJacobianColorvec
ArrayInterface.matrix_colors
— Functionmatrix_colors(A, alg::ColoringAlgorithm = GreedyD1Color();
partition_by_rows::Bool = false)
Return the colorvec vector for the matrix A using the chosen coloring algorithm. If a known analytical solution exists, that is used instead. The coloring defaults to a greedy distance-1 coloring.
Note that if A isa SparseMatrixCSC, the sparsity pattern is defined by structural nonzeroes, ie includes explicitly stored zeros.
If ArrayInterface.fast_matrix_colors(A)
is true, then uses ArrayInterface.matrix_colors(A)
to compute the matrix colors.
SparseDiffTools.JacVec
— FunctionJacVec(f, u, [p, t]; fu = nothing, autodiff = AutoForwardDiff(), tag = DeivVecTag(),
kwargs...)
Returns SciMLOperators.FunctionOperator which computes jacobian-vector product df/du * v
.
For non-square jacobians with inplace f
, fu
must be specified, else JacVec
assumes a square jacobian.
L = JacVec(f, u)
L * v # = df/du * v
mul!(w, L, v) # = df/du * v
L(v, p, t) # = df/dw * v
L(x, v, p, t) # = df/dw * v
Allowed Function Signatures for f
For Out of Place Functions:
f(u, p, t) # t !== nothing
f(u, p) # p !== nothing
f(u) # Otherwise
For In Place Functions:
f(du, u, p, t) # t !== nothing
f(du, u, p) # p !== nothing
f(du, u) # Otherwise
SparseDiffTools.VecJac
— FunctionVecJac(f, u, [p, t]; fu = nothing, autodiff = AutoFiniteDiff())
Returns SciMLOperators.FunctionOperator which computes vector-jacobian product (df/du)ᵀ * v
.
For non-square jacobians with inplace f
, fu
must be specified, else VecJac
assumes a square jacobian.
L = VecJac(f, u)
L * v # = (df/du)ᵀ * v
mul!(w, L, v) # = (df/du)ᵀ * v
L(v, p, t; VJP_input = w) # = (df/du)ᵀ * v
L(x, v, p, t; VJP_input = w) # = (df/du)ᵀ * v
Allowed Function Signatures for f
For Out of Place Functions:
f(u, p, t) # t !== nothing
f(u, p) # p !== nothing
f(u) # Otherwise
For In Place Functions:
f(du, u, p, t) # t !== nothing
f(du, u, p) # p !== nothing
f(du, u) # Otherwise
SparseDiffTools._cols_by_rows
— Method_cols_by_rows(rows_index,cols_index)
Returns a vector of rows where each row contains a vector of its column indices.
SparseDiffTools._rows_by_cols
— Method_rows_by_cols(rows_index,cols_index)
Returns a vector of columns where each column contains a vector of its row indices.
SparseDiffTools.color_graph
— Methodcolor_graph(g::Graphs.AbstractGraphs, ::AcyclicColoring)
Returns a coloring vector following the acyclic coloring rules (1) the coloring corresponds to a distance-1 coloring, and (2) vertices in every cycle of the graph are assigned at least three distinct colors. This variant of coloring is called acyclic since every subgraph induced by vertices assigned any two colors is a collection of trees—and hence is acyclic.
Reference: Gebremedhin AH, Manne F, Pothen A. New Acyclic and Star Coloring Algorithms with Application to Computing Hessians
SparseDiffTools.color_graph
— Methodcolor_graph(g::Graphs.AbstractGraph, ::BacktrackingColor)
Return a tight, distance-1 coloring of graph g using the minimum number of colors possible (i.e. the chromatic number of graph, χ(g)
)
SparseDiffTools.color_graph
— Methodcolor_graph(g::Graphs.AbstractGraph, ::GreedyStar1Color)
Find a coloring of a given input graph such that no two vertices connected by an edge have the same color using greedy approach. The number of colors used may be equal or greater than the chromatic number χ(G)
of the graph.
A star coloring is a special type of distance - 1 coloring, For a coloring to be called a star coloring, it must satisfy two conditions:
- every pair of adjacent vertices receives distinct colors
(a distance-1 coloring)
- For any vertex v, any color that leads to a two-colored path
involving v and three other vertices is impermissible for v. In other words, every path on four vertices uses at least three colors.
Reference: Gebremedhin AH, Manne F, Pothen A. What color is your Jacobian? Graph coloring for computing derivatives. SIAM review. 2005;47(4):629-705.
SparseDiffTools.color_graph
— Methodcolor_graph(g::Graphs.AbstractGraph, ::GreedyStar2Color)
Find a coloring of a given input graph such that no two vertices connected by an edge have the same color using greedy approach. The number of colors used may be equal or greater than the chromatic number χ(G)
of the graph.
A star coloring is a special type of distance - 1 coloring, For a coloring to be called a star coloring, it must satisfy two conditions:
- every pair of adjacent vertices receives distinct colors
(a distance-1 coloring)
- For any vertex v, any color that leads to a two-colored path
involving v and three other vertices is impermissible for v. In other words, every path on four vertices uses at least three colors.
Reference: Gebremedhin AH, Manne F, Pothen A. What color is your Jacobian? Graph coloring for computing derivatives. SIAM review. 2005;47(4):629-705.
TODO: add text explaining the difference between star1 and star2
SparseDiffTools.color_graph
— Methodcolor_graph(G::VSafeGraph,::ContractionColor)
Find a coloring of the graph g such that no two vertices connected by an edge have the same color.
SparseDiffTools.color_graph
— Methodcolor_graph(g::VSafeGraph, alg::GreedyD1Color)
Find a coloring of a given input graph such that no two vertices connected by an edge have the same color using greedy approach. The number of colors used may be equal or greater than the chromatic number χ(G)
of the graph.
SparseDiffTools.contract!
— Methodcontract!(g, y, x)
Contract the vertex y to x, both of which belong to graph G, that is delete vertex y and join x with the neighbors of y if they are not already connected with an edge.
SparseDiffTools.find
— Methodfind(w::Integer, x::Integer, g::Graphs.AbstractGraph,
two_colored_forest::DisjointSets{<:Integer})
Returns the root of the disjoint set to which the edge connecting vertices w and x in the graph g belongs to
SparseDiffTools.find_edge_index
— Methodfind_edge(g::Graphs.AbstractGraph, v::Integer, w::Integer)
Returns an integer equivalent to the index of the edge connecting the vertices v and w in the graph g
SparseDiffTools.free_colors
— Methodfree_colors(x::Integer,
A::AbstractVector{<:Integer},
colors::AbstractVector{<:Integer},
F::Array{Integer,1},
g::Graphs.AbstractGraph,
opt::Integer)
Returns set of free colors of x which are less than optimal chromatic number (opt)
Arguments:
x: Vertex who's set of free colors is to be calculated A: List of vertices of graph g sorted in non-increasing order of degree colors: colors[i] stores the number of distinct colors used in the coloring of vertices A[0], A[1]... A[i-1] F: F[i] stores the color of vertex i g: Graph to be colored opt: Current optimal number of colors to be used in the coloring of graph g
SparseDiffTools.grow_star!
— Methodgrow_star!(two_colored_forest::DisjointSets{<:Integer},
first_neighbor::AbstractVector{<:Tuple{Integer, Integer}}, v::Integer, w::Integer,
g::Graphs.AbstractGraph, color::AbstractVector{<:Integer})
Grow a 2-colored star after assigning a new color to the previously uncolored vertex v, by comparing it with the adjacent vertex w. Disjoint set is used to store stars in sets, which are identified through key edges present in g.
SparseDiffTools.init_jacobian
— Functioninit_jacobian(cache::AbstractMaybeSparseJacobianCache)
Initialize the Jacobian based on the cache. Uses sparse jacobians if possible.
This function doesn't alias the provided jacobian prototype. It always initializes a fresh jacobian that can be mutated without any side effects.
SparseDiffTools.init_jacobian
— Methodinit_jacobian(cache::AbstractMaybeSparseJacobianCache;
preserve_immutable::Val = Val(false))
Initialize the Jacobian based on the cache. Uses sparse jacobians if possible.
If preserve_immutable
is true
, then the Jacobian returned might be immutable, this is relevant if the inputs are immutable like StaticArrays
.
SparseDiffTools.insert_new_tree!
— Methodinsert_new_tree!(two_colored_forest::DisjointSets{<:Integer}, v::Integer,
w::Integer, g::Graphs.AbstractGraph
creates a new singleton set in the disjoint set 'twocoloredforest' consisting of the edge connecting v and w in the graph g
SparseDiffTools.least_index
— Methodleast_index(F::AbstractVector{<:Integer}, A::AbstractVector{<:Integer}, opt::Integer)
Returns least index i such that color of vertex A[i] is equal to opt
(optimal chromatic number)
SparseDiffTools.length_common_neighbor
— Methodlength_common_neighbor(g, z, x)
Find the number of vertices that share an edge with both the vertices z and x belonging to the graph g.
SparseDiffTools.matrix2graph
— Functionmatrix2graph(sparse_matrix, [partition_by_rows::Bool=true])
A utility function to generate a graph from input sparse matrix, columns are represented with vertices and 2 vertices are connected with an edge only if the two columns are mutually orthogonal.
Note that the sparsity pattern is defined by structural nonzeroes, ie includes explicitly stored zeros.
SparseDiffTools.max_degree_vertex
— Methodmax_degree_vertex(G, nn)
Find the vertex in the group nn of vertices belonging to the graph G which has the highest degree.
SparseDiffTools.max_degree_vertex
— Methodmax_degree_vertex(G)
Find the vertex in graph with highest degree.
SparseDiffTools.merge_trees!
— Methodmerge_trees!(two_colored_forest::DisjointSets{<:Integer}, v::Integer, w::Integer,
x::Integer, g::Graphs.AbstractGraph)
Subroutine to merge trees present in the disjoint set which have a common edge.
SparseDiffTools.min_index
— Methodmin_index(forbidden_colors::AbstractVector{<:Integer}, v::Integer)
Returns min{i > 0 such that forbidden_colors[i] != v}
SparseDiffTools.non_neighbors
— Methodnon_neighbors(G, x)
Find the set of vertices belonging to the graph G which do not share an edge with the vertex x.
SparseDiffTools.prevent_cycle!
— Methodprevent_cycle!(first_visit_to_tree::AbstractVector{<:Tuple{Integer, Integer}},
forbidden_colors::AbstractVector{<:Integer}, v::Integer, w::Integer, x::Integer,
g::Graphs.AbstractGraph, two_colored_forest::DisjointSets{<:Integer},
color::AbstractVector{<:Integer})
Subroutine to avoid generation of 2-colored cycle due to coloring of vertex v, which is adjacent to vertices w and x in graph g. Disjoint set is used to store the induced 2-colored subgraphs/trees where the id of set is an integer representing an edge of graph 'g'
SparseDiffTools.remove_higher_colors
— Methodremove_higher_colors(U::AbstractVector{<:Integer}, opt::Integer)
Remove all the colors which are greater than or equal to the opt
(optimal chromatic number) from the set of colors U
SparseDiffTools.sort_by_degree
— Methodsort_by_degree(g::Graphs.AbstractGraph)
Returns a list of the vertices of graph g sorted in non-increasing order of their degrees
SparseDiffTools.sparse_jacobian!
— Functionsparse_jacobian!(J::AbstractMatrix, ad, cache::AbstractMaybeSparseJacobianCache, f, x)
sparse_jacobian!(J::AbstractMatrix, ad, cache::AbstractMaybeSparseJacobianCache, f!, fx,
x)
Inplace update the matrix J
with the Jacobian of f
at x
using the AD backend ad
.
cache
is the cache object returned by sparse_jacobian_cache
.
SparseDiffTools.sparse_jacobian!
— Methodsparse_jacobian!(J::AbstractMatrix, ad::AbstractADType, sd::AbstractSparsityDetection,
f, x; fx=nothing)
sparse_jacobian!(J::AbstractMatrix, ad::AbstractADType, sd::AbstractSparsityDetection,
f!, fx, x)
Sequentially calls sparse_jacobian_cache
and sparse_jacobian!
to compute the Jacobian of f
at x
. Use this if the jacobian for f
is computed exactly once. In all other cases, use sparse_jacobian_cache
once to generate the cache and use sparse_jacobian!
with the same cache to compute the jacobian.
SparseDiffTools.sparse_jacobian
— Methodsparse_jacobian(ad::AbstractADType, cache::AbstractMaybeSparseJacobianCache, f, x)
sparse_jacobian(ad::AbstractADType, cache::AbstractMaybeSparseJacobianCache, f!, fx, x)
Use the sparsity detection cache
for computing the sparse Jacobian. This allocates a new Jacobian at every function call.
If x
is a StaticArray, then this function tries to use a non-allocating implementation for the jacobian computation. This is possible only for a limited backends currently.
SparseDiffTools.sparse_jacobian
— Methodsparse_jacobian(ad::AbstractADType, sd::AbstractMaybeSparsityDetection, f, x; fx=nothing)
sparse_jacobian(ad::AbstractADType, sd::AbstractMaybeSparsityDetection, f!, fx, x)
Sequentially calls sparse_jacobian_cache
and sparse_jacobian!
to compute the Jacobian of f
at x
. Use this if the jacobian for f
is computed exactly once. In all other cases, use sparse_jacobian_cache
once to generate the cache and use sparse_jacobian!
with the same cache to compute the jacobian.
If x
is a StaticArray, then this function tries to use a non-allocating implementation for the jacobian computation. This is possible only for a limited backends currently.
SparseDiffTools.sparse_jacobian_cache
— Methodsparse_jacobian_cache(ad::AbstractADType, sd::AbstractSparsityDetection, f, x;
fx=nothing)
sparse_jacobian_cache(ad::AbstractADType, sd::AbstractSparsityDetection, f!, fx, x)
Takes the underlying AD backend ad
, sparsity detection algorithm sd
, function f
, and input x
and returns a cache object that can be used to compute the Jacobian.
If fx
is not specified, it will be computed by calling f(x)
.
Returns
A cache for computing the Jacobian of type AbstractMaybeSparseJacobianCache
.
SparseDiffTools.uncolor_all!
— Methoduncolor_all(F::AbstractVector{<:Integer}, A::AbstractVector{<:Integer}, start::Integer)
Uncolors all vertices A[i] where i is greater than or equal to start
SparseDiffTools.uncolored_vertex_of_maximal_degree
— Methoduncolored_vertex_of_maximal_degree(A::AbstractVector{<:Integer},F::AbstractVector{<:Integer})
Returns an uncolored vertex from the partially colored graph which has the highest degree
SparseDiffTools.vertex_degree
— Methodvertex_degree(g, z)
Find the degree of the vertex z which belongs to the graph g.