Skip to content

Commit

Permalink
Merge 5263f5b into 320cd9c
Browse files Browse the repository at this point in the history
  • Loading branch information
ChrisRackauckas committed Apr 2, 2019
2 parents 320cd9c + 5263f5b commit 84a890c
Show file tree
Hide file tree
Showing 4 changed files with 47 additions and 0 deletions.
10 changes: 10 additions & 0 deletions src/initdt.jl
Expand Up @@ -102,6 +102,11 @@
f₁ .= ftmp
end

# Constant zone before callback
# Just return first guess
# Avoids AD issues
f₀ == f₁ && return 100dt₀

@. tmp = (f₁-f₀)/sk*oneunit_tType
d₂ = internalnorm(tmp,t)/dt₀*oneunit_tType
# Hairer has d₂ = sqrt(sum(abs2,tmp))/dt₀, note the lack of norm correction
Expand Down Expand Up @@ -143,6 +148,11 @@ end
u₁ = @. u0 + dt₀_tdir * f₀
f₁ = f(u₁,p,t+dt₀_tdir)

# Constant zone before callback
# Just return first guess
# Avoids AD issues
f₀ == f₁ && return 100dt₀

d₂ = internalnorm((f₁ .- f₀) ./ sk .* oneunit_tType,t) / dt₀ * oneunit_tType

max_d₁d₂ = max(d₁, d₂)
Expand Down
1 change: 1 addition & 0 deletions test/REQUIRE
Expand Up @@ -5,3 +5,4 @@ StaticArrays
DiffEqCallbacks
DiffEqOperators
SafeTestsets
Calculus
35 changes: 35 additions & 0 deletions test/ad_tests.jl
@@ -0,0 +1,35 @@
using ParameterizedFunctions, Test
using OrdinaryDiffEq, Calculus, ForwardDiff

f = @ode_def begin
dx = -a
dy = b
end a b

cb = ContinuousCallback((u,t,i) -> u[1], (integrator)->(println("Stopped.");integrator.p[2]=zero(integrator.p[2])))
function test_f(p)
prob = ODEProblem(f,eltype(p).([1.0,0.0]),eltype(p).((0.0,1.0)),copy(p))
integrator = init(prob,Tsit5(),abstol=1e-14,reltol=1e-14,callback=cb)
step!(integrator)
solve!(integrator).u[end]
end
p = [2.0, 1.0]
findiff = Calculus.finite_difference_jacobian(test_f,p)
fordiff = ForwardDiff.jacobian(test_f,p)
@test findiff fordiff

f2 = @ode_def begin
dx = -x
dy = b
end a b

function test_f2(p)
prob = ODEProblem(f2,eltype(p).([1.0,0.0]),eltype(p).((0.0,1.0)),copy(p))
integrator = init(prob,Tsit5(),abstol=1e-14,reltol=1e-14,callback=cb)
step!(integrator)
solve!(integrator).u[end]
end
p = [2.0, 1.0]
findiff = Calculus.finite_difference_jacobian(test_f2,p)
fordiff = ForwardDiff.jacobian(test_f2,p)
@test findiff fordiff
1 change: 1 addition & 0 deletions test/runtests.jl
Expand Up @@ -33,6 +33,7 @@ if group == "All" || group == "Interface"
@time @safetestset "Derivative Utilities Tests" begin include("utility_tests.jl") end
@time @safetestset "Discrete Callback Dual Tests" begin include("discrete_callback_dual_test.jl") end
@time @safetestset "DEStats Tests" begin include("destats_tests.jl") end
@time @safetestset "AD Tests" begin include("ad_tests.jl") end
end

if group == "All" || group == "Integrators"
Expand Down

0 comments on commit 84a890c

Please sign in to comment.