Skip to content

Commit d66a379

Browse files
committed
Add support for StaticArrays >= 1.7 (#703)
1 parent 0feb64c commit d66a379

File tree

5 files changed

+78
-18
lines changed

5 files changed

+78
-18
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -26,7 +26,7 @@ LogExpFunctions = "0.3"
2626
NaNMath = "0.2.2, 0.3, 1"
2727
Preferences = "1"
2828
SpecialFunctions = "0.8, 0.9, 0.10, 1.0, 2"
29-
StaticArrays = "0.8.3, 0.9, 0.10, 0.11, 0.12, ~1.0, ~1.1, ~1.2, ~1.3, ~1.4, ~1.5, ~1.6"
29+
StaticArrays = "0.8.3, 0.9, 0.10, 0.11, 0.12, 1.0"
3030
julia = "1"
3131

3232
[extensions]

ext/ForwardDiffStaticArraysExt.jl

Lines changed: 7 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -7,7 +7,7 @@ using ForwardDiff: Dual, partials, GradientConfig, JacobianConfig, HessianConfig
77
gradient, hessian, jacobian, gradient!, hessian!, jacobian!,
88
extract_gradient!, extract_jacobian!, extract_value!,
99
vector_mode_gradient, vector_mode_gradient!,
10-
vector_mode_jacobian, vector_mode_jacobian!, valtype, value, _lyap_div!
10+
vector_mode_jacobian, vector_mode_jacobian!, valtype, value
1111
using DiffResults: DiffResult, ImmutableDiffResult, MutableDiffResult
1212

1313
@generated function dualize(::Type{T}, x::StaticArray) where T
@@ -23,19 +23,17 @@ end
2323

2424
@inline static_dual_eval(::Type{T}, f, x::StaticArray) where T = f(dualize(T, x))
2525

26+
# To fix method ambiguity issues:
2627
function LinearAlgebra.eigvals(A::Symmetric{<:Dual{Tg,T,N}, <:StaticArrays.StaticMatrix}) where {Tg,T<:Real,N}
27-
λ,Q = eigen(Symmetric(value.(parent(A))))
28-
parts = ntuple(j -> diag(Q' * getindex.(partials.(A), j) * Q), N)
29-
Dual{Tg}.(λ, tuple.(parts...))
28+
return ForwardDiff._eigvals(A)
3029
end
31-
3230
function LinearAlgebra.eigen(A::Symmetric{<:Dual{Tg,T,N}, <:StaticArrays.StaticMatrix}) where {Tg,T<:Real,N}
33-
λ = eigvals(A)
34-
_,Q = eigen(Symmetric(value.(parent(A))))
35-
parts = ntuple(j -> Q*ForwardDiff._lyap_div!(Q' * getindex.(partials.(A), j) * Q - Diagonal(getindex.(partials.(λ), j)), value.(λ)), N)
36-
Eigen(λ,Dual{Tg}.(Q, tuple.(parts...)))
31+
return ForwardDiff._eigen(A)
3732
end
3833

34+
# For `MMatrix` we can use the in-place method
35+
ForwardDiff._lyap_div!!(A::StaticArrays.MMatrix, λ::AbstractVector) = ForwardDiff._lyap_div!(A, λ)
36+
3937
# Gradient
4038
@inline ForwardDiff.gradient(f::F, x::StaticArray) where F = vector_mode_gradient(f, x)
4139
@inline ForwardDiff.gradient(f::F, x::StaticArray, cfg::GradientConfig) where F = gradient(f, x)

src/dual.jl

Lines changed: 25 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -702,7 +702,11 @@ end
702702
# Symmetric eigvals #
703703
#-------------------#
704704

705-
function LinearAlgebra.eigvals(A::Symmetric{<:Dual{Tg,T,N}}) where {Tg,T<:Real,N}
705+
# To be able to reuse this default definition in the StaticArrays extension
706+
# (has to be re-defined to avoid method ambiguity issues)
707+
# we forward the call to an internal method that can be shared and reused
708+
LinearAlgebra.eigvals(A::Symmetric{<:Dual{Tg,T,N}}) where {Tg,T<:Real,N} = _eigvals(A)
709+
function _eigvals(A::Symmetric{<:Dual{Tg,T,N}}) where {Tg,T<:Real,N}
706710
λ,Q = eigen(Symmetric(value.(parent(A))))
707711
parts = ntuple(j -> diag(Q' * getindex.(partials.(A), j) * Q), N)
708712
Dual{Tg}.(λ, tuple.(parts...))
@@ -720,8 +724,19 @@ function LinearAlgebra.eigvals(A::SymTridiagonal{<:Dual{Tg,T,N}}) where {Tg,T<:R
720724
Dual{Tg}.(λ, tuple.(parts...))
721725
end
722726

723-
# A ./ (λ - λ') but with diag special cased
724-
function _lyap_div!(A, λ)
727+
# A ./ (λ' .- λ) but with diag special cased
728+
# Default out-of-place method
729+
function _lyap_div!!(A::AbstractMatrix, λ::AbstractVector)
730+
return map(
731+
(a, b, idx) -> a / (idx[1] == idx[2] ? oneunit(b) : b),
732+
A,
733+
λ' .- λ,
734+
CartesianIndices(A),
735+
)
736+
end
737+
# For `Matrix` (and e.g. `StaticArrays.MMatrix`) we can use an in-place method
738+
_lyap_div!!(A::Matrix, λ::AbstractVector) = _lyap_div!(A, λ)
739+
function _lyap_div!(A::AbstractMatrix, λ::AbstractVector)
725740
for (j,μ) in enumerate(λ), (k,λ) in enumerate(λ)
726741
if k j
727742
A[k,j] /= μ - λ
@@ -730,17 +745,21 @@ function _lyap_div!(A, λ)
730745
A
731746
end
732747

733-
function LinearAlgebra.eigen(A::Symmetric{<:Dual{Tg,T,N}}) where {Tg,T<:Real,N}
748+
# To be able to reuse this default definition in the StaticArrays extension
749+
# (has to be re-defined to avoid method ambiguity issues)
750+
# we forward the call to an internal method that can be shared and reused
751+
LinearAlgebra.eigen(A::Symmetric{<:Dual{Tg,T,N}}) where {Tg,T<:Real,N} = _eigen(A)
752+
function _eigen(A::Symmetric{<:Dual{Tg,T,N}}) where {Tg,T<:Real,N}
734753
λ = eigvals(A)
735754
_,Q = eigen(Symmetric(value.(parent(A))))
736-
parts = ntuple(j -> Q*_lyap_div!(Q' * getindex.(partials.(A), j) * Q - Diagonal(getindex.(partials.(λ), j)), value.(λ)), N)
755+
parts = ntuple(j -> Q*_lyap_div!!(Q' * getindex.(partials.(A), j) * Q - Diagonal(getindex.(partials.(λ), j)), value.(λ)), N)
737756
Eigen(λ,Dual{Tg}.(Q, tuple.(parts...)))
738757
end
739758

740759
function LinearAlgebra.eigen(A::SymTridiagonal{<:Dual{Tg,T,N}}) where {Tg,T<:Real,N}
741760
λ = eigvals(A)
742761
_,Q = eigen(SymTridiagonal(value.(parent(A))))
743-
parts = ntuple(j -> Q*_lyap_div!(Q' * getindex.(partials.(A), j) * Q - Diagonal(getindex.(partials.(λ), j)), value.(λ)), N)
762+
parts = ntuple(j -> Q*_lyap_div!!(Q' * getindex.(partials.(A), j) * Q - Diagonal(getindex.(partials.(λ), j)), value.(λ)), N)
744763
Eigen(λ,Dual{Tg}.(Q, tuple.(parts...)))
745764
end
746765

test/JacobianTest.jl

Lines changed: 6 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -239,8 +239,12 @@ end
239239
@test ForwardDiff.jacobian(ev1, x0) Calculus.finite_difference_jacobian(ev1, x0)
240240
ev2(x) = eigen(SymTridiagonal(x, x[1:end-1])).vectors[:,1]
241241
@test ForwardDiff.jacobian(ev2, x0) Calculus.finite_difference_jacobian(ev2, x0)
242-
x0_static = SVector{2}(x0)
243-
@test ForwardDiff.jacobian(ev1, x0_static) Calculus.finite_difference_jacobian(ev1, x0)
242+
243+
x0_svector = SVector{2}(x0)
244+
@test ForwardDiff.jacobian(ev1, x0_svector) Calculus.finite_difference_jacobian(ev1, x0)
245+
246+
x0_mvector = MVector{2}(x0)
247+
@test ForwardDiff.jacobian(ev1, x0_mvector) Calculus.finite_difference_jacobian(ev1, x0)
244248
end
245249

246250
@testset "type stability" begin

test/MiscTest.jl

Lines changed: 39 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -6,7 +6,9 @@ using Test
66
using ForwardDiff
77
using DiffTests
88
using SparseArrays: sparse
9+
using StaticArrays
910
using IrrationalConstants
11+
using LinearAlgebra
1012

1113
include(joinpath(dirname(@__FILE__), "utils.jl"))
1214

@@ -177,4 +179,41 @@ end
177179
@test ForwardDiff.derivative(x -> rem2pi(x, RoundUp), rand()) == 1
178180
@test ForwardDiff.derivative(x -> rem2pi(x, RoundDown), rand()) == 1
179181

182+
@testset "_lyap_div!!" begin
183+
# In-place version for `Matrix`
184+
A = rand(3, 3)
185+
Acopy = copy(A)
186+
λ = rand(3)
187+
B = @inferred(ForwardDiff._lyap_div!!(A, λ))
188+
@test B === A
189+
@test B[diagind(B)] == Acopy[diagind(Acopy)]
190+
no_diag(X) = [X[i] for i in eachindex(X) if !(i in diagind(X))]
191+
@test no_diag(B) == no_diag(Acopy ./' .- λ))
192+
193+
# Immutable static arrays
194+
A_smatrix = SMatrix{3,3}(Acopy)
195+
λ_svector = SVector{3}(λ)
196+
B_smatrix = @inferred(ForwardDiff._lyap_div!!(A_smatrix, λ_svector))
197+
@test B_smatrix !== A_smatrix
198+
@test B_smatrix isa SMatrix{3,3}
199+
@test B_smatrix == B
200+
λ_mvector = MVector{3}(λ)
201+
B_smatrix = @inferred(ForwardDiff._lyap_div!!(A_smatrix, λ_mvector))
202+
@test B_smatrix !== A_smatrix
203+
@test B_smatrix isa SMatrix{3,3}
204+
@test B_smatrix == B
205+
206+
# Mutable static arrays
207+
A_mmatrix = MMatrix{3,3}(Acopy)
208+
λ_svector = SVector{3}(λ)
209+
B_mmatrix = @inferred(ForwardDiff._lyap_div!!(A_mmatrix, λ_svector))
210+
@test B_mmatrix === A_mmatrix
211+
@test B_mmatrix == B
212+
A_mmatrix = MMatrix{3,3}(Acopy)
213+
λ_mvector = MVector{3}(λ)
214+
B_mmatrix = @inferred(ForwardDiff._lyap_div!!(A_mmatrix, λ_mvector))
215+
@test B_mmatrix === A_mmatrix
216+
@test B_mmatrix == B
217+
end
218+
180219
end # module

0 commit comments

Comments
 (0)