Skip to content

ForwardDiff.gradient and ForwardDiff.jacobian don't handle functions with ForwardDiff.Dual output correctly #769

@devmotion

Description

@devmotion

MWE:

julia> using ForwardDiff, StaticArrays

julia> struct TestTag end

julia> struct OuterTestTag end

julia> ForwardDiff.:(::Type{TestTag}, ::Type{OuterTestTag}) = true

julia> ForwardDiff.:(::Type{OuterTestTag}, ::Type{<:ForwardDiff.Tag}) = true

julia> x = [ForwardDiff.Dual{OuterTestTag}(ForwardDiff.Dual{TestTag}(1.3, 2.1), ForwardDiff.Dual{TestTag}(0.3, -2.4))];

julia> ForwardDiff.derivative(ForwardDiff.value, only(x))
Dual{OuterTestTag}(Dual{TestTag}(0.0,0.0),Dual{TestTag}(0.0,0.0))

julia> ForwardDiff.gradient(x -> ForwardDiff.value(only(x)), x)
ERROR: MethodError: no method matching Float64(::ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{TestTag, Float64, 1}, 1})
The type `Float64` exists, but no method is defined for this combination of argument types when trying to construct it.

Closest candidates are:
  (::Type{T})(::Real, ::RoundingMode) where T<:AbstractFloat
   @ Base rounding.jl:265
  (::Type{T})(::T) where T<:Number
   @ Core boot.jl:900
  Float64(::IrrationalConstants.Inv2π)
   @ IrrationalConstants ~/.julia/packages/IrrationalConstants/lWTip/src/macro.jl:131
  ...

Stacktrace:
 [1] convert
   @ ~/.julia/packages/ForwardDiff/Or6Qh/src/dual.jl:467 [inlined]
 [2] setindex!(A::Vector{ForwardDiff.Dual{…}}, x::ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{…}, 1}, i::Int64)
   @ Base ./array.jl:987
 [3] extract_gradient!
   @ ~/.julia/packages/ForwardDiff/Or6Qh/src/gradient.jl:69 [inlined]
 [4] vector_mode_gradient
   @ ~/.julia/packages/ForwardDiff/Or6Qh/src/gradient.jl:101 [inlined]
 [5] gradient
   @ ~/.julia/packages/ForwardDiff/Or6Qh/src/gradient.jl:20 [inlined]
 [6] gradient(f::var"#9#10", x::Vector{…}, cfg::ForwardDiff.GradientConfig{…})
   @ ForwardDiff ~/.julia/packages/ForwardDiff/Or6Qh/src/gradient.jl:17
 [7] gradient(f::var"#9#10", x::Vector{ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{TestTag, Float64, 1}, 1}})
   @ ForwardDiff ~/.julia/packages/ForwardDiff/Or6Qh/src/gradient.jl:17
 [8] top-level scope
   @ REPL[11]:1
Some type information was truncated. Use `show(err)` to see complete types.

julia> ForwardDiff.gradient(x -> ForwardDiff.value(only(x)), SVector{1}(x))
ERROR: MethodError: no method matching Float64(::ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{TestTag, Float64, 1}, 1})
The type `Float64` exists, but no method is defined for this combination of argument types when trying to construct it.

Closest candidates are:
  (::Type{T})(::Real, ::RoundingMode) where T<:AbstractFloat
   @ Base rounding.jl:265
  (::Type{T})(::T) where T<:Number
   @ Core boot.jl:900
  Float64(::IrrationalConstants.Inv2π)
   @ IrrationalConstants ~/.julia/packages/IrrationalConstants/lWTip/src/macro.jl:131
  ...

Stacktrace:
 [1] convert
   @ ~/.julia/packages/ForwardDiff/Or6Qh/src/dual.jl:467 [inlined]
 [2] macro expansion
   @ ~/.julia/packages/StaticArraysCore/7xxEJ/src/StaticArraysCore.jl:88 [inlined]
 [3] convert_ntuple
   @ ~/.julia/packages/StaticArraysCore/7xxEJ/src/StaticArraysCore.jl:84 [inlined]
 [4] SVector{1, ForwardDiff.Dual{TestTag, Float64, 1}}(x::Tuple{ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{…}, 1}})
   @ StaticArraysCore ~/.julia/packages/StaticArraysCore/7xxEJ/src/StaticArraysCore.jl:120
 [5] macro expansion
   @ ~/.julia/packages/ForwardDiff/Or6Qh/ext/ForwardDiffStaticArraysExt.jl:51 [inlined]
 [6] extract_gradient
   @ ~/.julia/packages/ForwardDiff/Or6Qh/ext/ForwardDiffStaticArraysExt.jl:46 [inlined]
 [7] vector_mode_gradient
   @ ~/.julia/packages/ForwardDiff/Or6Qh/ext/ForwardDiffStaticArraysExt.jl:57 [inlined]
 [8] gradient(f::var"#13#14", x::SVector{1, ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{TestTag, Float64, 1}, 1}})
   @ ForwardDiffStaticArraysExt ~/.julia/packages/ForwardDiff/Or6Qh/ext/ForwardDiffStaticArraysExt.jl:38
 [9] top-level scope
   @ REPL[12]:1
Some type information was truncated. Use `show(err)` to see complete types.

julia> ForwardDiff.jacobian(x -> map(ForwardDiff.value, x), x)
ERROR: MethodError: no method matching Float64(::ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{TestTag, Float64, 1}, 1})
The type `Float64` exists, but no method is defined for this combination of argument types when trying to construct it.

Closest candidates are:
  (::Type{T})(::Real, ::RoundingMode) where T<:AbstractFloat
   @ Base rounding.jl:265
  (::Type{T})(::T) where T<:Number
   @ Core boot.jl:900
  Float64(::IrrationalConstants.Inv2π)
   @ IrrationalConstants ~/.julia/packages/IrrationalConstants/lWTip/src/macro.jl:131
  ...

Stacktrace:
  [1] convert
    @ ~/.julia/packages/ForwardDiff/Or6Qh/src/dual.jl:467 [inlined]
  [2] setindex!
    @ ./array.jl:994 [inlined]
  [3] setindex!
    @ ./multidimensional.jl:704 [inlined]
  [4] macro expansion
    @ ./broadcast.jl:973 [inlined]
  [5] macro expansion
    @ ./simdloop.jl:77 [inlined]
  [6] copyto!
    @ ./broadcast.jl:972 [inlined]
  [7] copyto!
    @ ./broadcast.jl:925 [inlined]
  [8] materialize!
    @ ./broadcast.jl:883 [inlined]
  [9] materialize!
    @ ./broadcast.jl:880 [inlined]
 [10] extract_jacobian!(::Type{…}, result::Matrix{…}, ydual::Vector{…}, n::Int64)
    @ ForwardDiff ~/.julia/packages/ForwardDiff/Or6Qh/src/jacobian.jl:100
 [11] vector_mode_jacobian(f::var"#15#16", x::Vector{…}, cfg::ForwardDiff.JacobianConfig{…})
    @ ForwardDiff ~/.julia/packages/ForwardDiff/Or6Qh/src/jacobian.jl:132
 [12] jacobian
    @ ~/.julia/packages/ForwardDiff/Or6Qh/src/jacobian.jl:22 [inlined]
 [13] jacobian(f::var"#15#16", x::Vector{…}, cfg::ForwardDiff.JacobianConfig{…})
    @ ForwardDiff ~/.julia/packages/ForwardDiff/Or6Qh/src/jacobian.jl:19
 [14] jacobian(f::var"#15#16", x::Vector{ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{TestTag, Float64, 1}, 1}})
    @ ForwardDiff ~/.julia/packages/ForwardDiff/Or6Qh/src/jacobian.jl:19
 [15] top-level scope
    @ REPL[13]:1
Some type information was truncated. Use `show(err)` to see complete types.

julia> ForwardDiff.jacobian(x -> map(ForwardDiff.value, x), SVector{1}(x))
ERROR: MethodError: no method matching Float64(::ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{TestTag, Float64, 1}, 1})
The type `Float64` exists, but no method is defined for this combination of argument types when trying to construct it.

Closest candidates are:
  (::Type{T})(::Real, ::RoundingMode) where T<:AbstractFloat
   @ Base rounding.jl:265
  (::Type{T})(::T) where T<:Number
   @ Core boot.jl:900
  Float64(::IrrationalConstants.Inv2π)
   @ IrrationalConstants ~/.julia/packages/IrrationalConstants/lWTip/src/macro.jl:131
  ...

Stacktrace:
 [1] convert
   @ ~/.julia/packages/ForwardDiff/Or6Qh/src/dual.jl:467 [inlined]
 [2] macro expansion
   @ ~/.julia/packages/StaticArraysCore/7xxEJ/src/StaticArraysCore.jl:88 [inlined]
 [3] convert_ntuple
   @ ~/.julia/packages/StaticArraysCore/7xxEJ/src/StaticArraysCore.jl:84 [inlined]
 [4] SMatrix{1, 1, ForwardDiff.Dual{…}, 1}(x::Tuple{ForwardDiff.Dual{…}})
   @ StaticArraysCore ~/.julia/packages/StaticArraysCore/7xxEJ/src/StaticArraysCore.jl:120
 [5] macro expansion
   @ ~/.julia/packages/ForwardDiff/Or6Qh/ext/ForwardDiffStaticArraysExt.jl:80 [inlined]
 [6] extract_jacobian
   @ ~/.julia/packages/ForwardDiff/Or6Qh/ext/ForwardDiffStaticArraysExt.jl:74 [inlined]
 [7] vector_mode_jacobian
   @ ~/.julia/packages/ForwardDiff/Or6Qh/ext/ForwardDiffStaticArraysExt.jl:86 [inlined]
 [8] jacobian(f::var"#19#20", x::SVector{1, ForwardDiff.Dual{OuterTestTag, ForwardDiff.Dual{TestTag, Float64, 1}, 1}})
   @ ForwardDiffStaticArraysExt ~/.julia/packages/ForwardDiff/Or6Qh/ext/ForwardDiffStaticArraysExt.jl:66
 [9] top-level scope
   @ REPL[14]:1
Some type information was truncated. Use `show(err)` to see complete types.

The problem is that it's incorrect to use valtype in

result = similar(x, valtype(ydual))
,
$(chunk_mode_gradient_expr(:(result = similar(x, valtype(ydual)))))
,
V = StaticArrays.similar_type(S, valtype($y))
,
V = StaticArrays.similar_type(S, valtype(eltype($ydual)), Size($M, $N))
,
result = similar(ydual, valtype(eltype(ydual)), length(ydual), length(x))
,
result = similar(ydual, valtype(eltype(ydual)), length(ydual), N)
, and
:(result = similar(ydual, valtype(eltype(ydual)), length(ydual), xlen)),
- in general, this does not return the type of ForwardDiff.value(T, dual) (the primal values of interest in these differentiations).

Metadata

Metadata

Assignees

No one assigned

    Labels

    Type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions