diff --git a/Project.toml b/Project.toml index 9eb2914c..2712d857 100644 --- a/Project.toml +++ b/Project.toml @@ -1,6 +1,6 @@ name = "ForwardDiff" uuid = "f6369f11-7733-5829-9624-2563aa707210" -version = "1.0.1" +version = "1.0.2" [deps] CommonSubexpressions = "bbf7d656-a473-5ed7-a52c-81e309532950" diff --git a/src/apiutils.jl b/src/apiutils.jl index 39608d6c..f401a3fc 100644 --- a/src/apiutils.jl +++ b/src/apiutils.jl @@ -72,16 +72,36 @@ end function seed!(duals::AbstractArray{Dual{T,V,N}}, x, seed::Partials{N,V} = zero(Partials{N,V})) where {T,V,N} - for idx in structural_eachindex(duals, x) - duals[idx] = Dual{T,V,N}(x[idx], seed) + if isbitstype(V) + for idx in structural_eachindex(duals, x) + duals[idx] = Dual{T,V,N}(x[idx], seed) + end + else + for idx in structural_eachindex(duals, x) + if isassigned(x, idx) + duals[idx] = Dual{T,V,N}(x[idx], seed) + else + Base._unsetindex!(duals, idx) + end + end end return duals end function seed!(duals::AbstractArray{Dual{T,V,N}}, x, seeds::NTuple{N,Partials{N,V}}) where {T,V,N} - for (i, idx) in zip(1:N, structural_eachindex(duals, x)) - duals[idx] = Dual{T,V,N}(x[idx], seeds[i]) + if isbitstype(V) + for (i, idx) in zip(1:N, structural_eachindex(duals, x)) + duals[idx] = Dual{T,V,N}(x[idx], seeds[i]) + end + else + for (i, idx) in zip(1:N, structural_eachindex(duals, x)) + if isassigned(x, idx) + duals[idx] = Dual{T,V,N}(x[idx], seeds[i]) + else + Base._unsetindex!(duals, idx) + end + end end return duals end @@ -90,8 +110,18 @@ function seed!(duals::AbstractArray{Dual{T,V,N}}, x, index, seed::Partials{N,V} = zero(Partials{N,V})) where {T,V,N} offset = index - 1 idxs = Iterators.drop(structural_eachindex(duals, x), offset) - for idx in idxs - duals[idx] = Dual{T,V,N}(x[idx], seed) + if isbitstype(V) + for idx in idxs + duals[idx] = Dual{T,V,N}(x[idx], seed) + end + else + for idx in idxs + if isassigned(x, idx) + duals[idx] = Dual{T,V,N}(x[idx], seed) + else + Base._unsetindex!(duals, idx) + end + end end return duals end @@ -100,8 +130,18 @@ function seed!(duals::AbstractArray{Dual{T,V,N}}, x, index, seeds::NTuple{N,Partials{N,V}}, chunksize = N) where {T,V,N} offset = index - 1 idxs = Iterators.drop(structural_eachindex(duals, x), offset) - for (i, idx) in zip(1:chunksize, idxs) - duals[idx] = Dual{T,V,N}(x[idx], seeds[i]) + if isbitstype(V) + for (i, idx) in zip(1:chunksize, idxs) + duals[idx] = Dual{T,V,N}(x[idx], seeds[i]) + end + else + for (i, idx) in zip(1:chunksize, idxs) + if isassigned(x, idx) + duals[idx] = Dual{T,V,N}(x[idx], seeds[i]) + else + Base._unsetindex!(duals, idx) + end + end end return duals end diff --git a/test/JacobianTest.jl b/test/JacobianTest.jl index 1e52f7fa..b010078f 100644 --- a/test/JacobianTest.jl +++ b/test/JacobianTest.jl @@ -279,4 +279,33 @@ end end end +# issues #436, #740 +@testset "BigFloat" begin + # Unassigned entries in the output + x = BigFloat.(1:9) + for chunksize in (1, 2, 9) + y = similar(x) + @test all(i -> !isassigned(y, i), eachindex(y)) + cfg = ForwardDiff.JacobianConfig(copyto!, y, x, ForwardDiff.Chunk{chunksize}()) + res = ForwardDiff.jacobian(copyto!, y, x, cfg) + @test y == x + @test res isa Matrix{BigFloat} + @test res == I + end + + # Unassigned (but unused) entry in the input and unassigned entries in the output + resize!(x, 10) + f = (y, x) -> copyto!(y, 1, x, 1, 9) + for chunksize in (1, 2, 10) + y = similar(x, 9) + @test all(i -> !isassigned(y, i), eachindex(y)) + cfg = ForwardDiff.JacobianConfig(f, y, x, ForwardDiff.Chunk{chunksize}()) + res = ForwardDiff.jacobian(f, y, x, cfg) + @test y == x[1:(end-1)] + @test res isa Matrix{BigFloat} + @test res[:, 1:(end-1)] == I + @test all(iszero, res[:, end]) + end +end + end # module