Structure and Sparsity Detection
Using the tracing system provided by Symbolics.jl expressions, Symbolics.jl can automatically detect the sparsity patterns of Julia functions efficiently. This functionality is described in more detail in the paper:
@article{gowda2019sparsity, title={Sparsity Programming: Automated Sparsity-Aware Optimizations in Differentiable Programming}, author={Gowda, Shashi and Ma, Yingbo and Churavy, Valentin and Edelman, Alan and Rackauckas, Christopher}, year={2019} }
Please cite this work if the functionality is used.
Sparsity Detection
#
Symbolics.jacobian_sparsity
— Function
jacobian_sparsity(
exprs::AbstractArray,
vars::AbstractArray
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}
Return the sparsity pattern of the Jacobian of an array of expressions with respect to an array of variable expressions.
Arguments
-
exprs
: an array of symbolic expressions. -
vars
: an array of symbolic variables.
Examples
julia> using Symbolics
julia> vars = @variables x₁ x₂;
julia> exprs = [2x₁, 3x₂, 4x₁ * x₂];
julia> Symbolics.jacobian_sparsity(exprs, vars)
3×2 SparseArrays.SparseMatrixCSC{Bool, Int64} with 4 stored entries:
1 ⋅
⋅ 1
1 1
jacobian_sparsity(
f!::Function,
output::AbstractArray,
input::AbstractArray,
args...;
kwargs...
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}
Return the sparsity pattern of the Jacobian of the mutating function f!
.
Arguments
-
f!
: an in-place functionf!(output, input, args...; kwargs...)
. -
output
: output array. -
input
: input array.
Examples
julia> using Symbolics
julia> f!(y, x) = y .= [x[2], 2x[1], 3x[1] * x[2]];
julia> output = Vector{Float64}(undef, 3);
julia> input = Vector{Float64}(undef, 2);
julia> Symbolics.jacobian_sparsity(f!, output, input)
3×2 SparseArrays.SparseMatrixCSC{Bool, Int64} with 4 stored entries:
⋅ 1
1 ⋅
1 1
#
Symbolics.hessian_sparsity
— Function
hessian_sparsity(
expr,
vars::AbstractVector;
full
) -> Union{SparseArrays.SparseMatrixCSC{Bool, Int64}, SparseArrays.SparseMatrixCSC{Int64, Int64}}
Return the sparsity pattern of the Hessian of an expression with respect to an array of variable expressions.
Arguments
-
expr
: a symbolic expression. -
vars
: a vector of symbolic variables.
Examples
julia> using Symbolics
julia> vars = @variables x₁ x₂;
julia> expr = 3x₁^2 + 4x₁ * x₂;
julia> Symbolics.hessian_sparsity(expr, vars)
2×2 SparseArrays.SparseMatrixCSC{Bool, Int64} with 3 stored entries:
1 1
1 ⋅
hessian_sparsity(
f::Function,
input::AbstractVector,
args...;
full,
kwargs...
) -> Union{SparseArrays.SparseMatrixCSC{Bool, Int64}, SparseArrays.SparseMatrixCSC{Int64, Int64}}
Return the sparsity pattern of the Hessian of the given function f
.
Arguments
Examples
julia> using Symbolics
julia> f(x) = 4x[1] * x[2] - 5x[2]^2;
julia> input = Vector{Float64}(undef, 2);
julia> Symbolics.hessian_sparsity(f, input)
2×2 SparseArrays.SparseMatrixCSC{Bool, Int64} with 3 stored entries:
⋅ 1
1 1