Skip to content

Remove inline and specialize the function evaluations #167

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 2 commits into from
Dec 8, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
14 changes: 7 additions & 7 deletions src/differentiation/compute_jacobian_ad.jl
Original file line number Diff line number Diff line change
Expand Up @@ -79,13 +79,13 @@ function generate_chunked_partials(x,colorvec,::Val{chunksize}) where chunksize
chunked_partials
end

@inline function forwarddiff_color_jacobian(f,
function forwarddiff_color_jacobian(f::F,
x::AbstractArray{<:Number};
colorvec = 1:length(x),
sparsity = nothing,
jac_prototype = nothing,
chunksize = nothing,
dx = sparsity === nothing && jac_prototype === nothing ? nothing : copy(x)) #if dx is nothing, we will estimate dx at the cost of a function call
dx = sparsity === nothing && jac_prototype === nothing ? nothing : copy(x)) where {F} #if dx is nothing, we will estimate dx at the cost of a function call

if sparsity === nothing && jac_prototype === nothing
cfg = if chunksize === nothing
Expand All @@ -95,7 +95,7 @@ end
ForwardDiff.JacobianConfig(f, x)
end
else
ForwardDiff.JacobianConfig(f, x, ForwardDiff.Chunk(getsize(chunksize)))
ForwardDiff.JacobianConfig(f, x, ForwardDiff.Chunk{getsize(chunksize)}())
end
return ForwardDiff.jacobian(f, x, cfg)
end
Expand All @@ -105,21 +105,21 @@ end
return forwarddiff_color_jacobian(f,x,ForwardColorJacCache(f,x,chunksize,dx=dx,colorvec=colorvec,sparsity=sparsity),jac_prototype)
end

@inline function forwarddiff_color_jacobian(J::AbstractArray{<:Number}, f,
function forwarddiff_color_jacobian(J::AbstractArray{<:Number}, f::F,
x::AbstractArray{<:Number};
colorvec = 1:length(x),
sparsity = nothing,
jac_prototype = nothing,
chunksize = nothing,
dx = similar(x, size(J, 1))) #dx kwarg can be used to avoid re-allocating dx every time
dx = similar(x, size(J, 1))) where {F} #dx kwarg can be used to avoid re-allocating dx every time
if sparsity === nothing && jac_prototype === nothing
cfg = chunksize === nothing ? ForwardDiff.JacobianConfig(f, x) : ForwardDiff.JacobianConfig(f, x, ForwardDiff.Chunk(getsize(chunksize)))
return ForwardDiff.jacobian(f, x, cfg)
end
return forwarddiff_color_jacobian(J,f,x,ForwardColorJacCache(f,x,chunksize,dx=dx,colorvec=colorvec,sparsity=sparsity))
end

function forwarddiff_color_jacobian(f,x::AbstractArray{<:Number},jac_cache::ForwardColorJacCache,jac_prototype=nothing)
function forwarddiff_color_jacobian(f::F,x::AbstractArray{<:Number},jac_cache::ForwardColorJacCache,jac_prototype=nothing) where F

if jac_prototype isa Nothing ? ArrayInterface.ismutable(x) : ArrayInterface.ismutable(jac_prototype)
# Whenever J is mutable, we mutate it to avoid allocations
Expand All @@ -136,7 +136,7 @@ function forwarddiff_color_jacobian(f,x::AbstractArray{<:Number},jac_cache::Forw
end

# When J is mutable, this version of forwarddiff_color_jacobian will mutate J to avoid allocations
function forwarddiff_color_jacobian(J::AbstractMatrix{<:Number},f,x::AbstractArray{<:Number},jac_cache::ForwardColorJacCache)
function forwarddiff_color_jacobian(J::AbstractMatrix{<:Number},f::F,x::AbstractArray{<:Number},jac_cache::ForwardColorJacCache) where F
t = jac_cache.t
dx = jac_cache.dx
p = jac_cache.p
Expand Down
2 changes: 1 addition & 1 deletion test/test_acyclic.jl
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,7 @@ using Graphs: SimpleGraph
using Test

using Random
Random.seed!(100)
Random.seed!(90)

# println("Starting acyclic coloring test...")
#= Test data =#
Expand Down