Skip to content

Commit 40aaab2

Browse files
Fix remake of ForwardDiffSensitivity
Fixes #1137
1 parent 16aad36 commit 40aaab2

File tree

2 files changed

+20
-2
lines changed

2 files changed

+20
-2
lines changed

src/forward_sensitivity.jl

Lines changed: 8 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -675,8 +675,14 @@ function SciMLBase.remake(
675675
{uType, tType, isinplace, P, F, K}
676676
_p = p === nothing ? parameter_values(prob) : p
677677
_f = f === nothing ? prob.f.f : f
678-
_u0 = u0 === nothing ? state_values(prob, 1:(prob.f.numindvar)) :
679-
u0[1:(prob.f.numindvar)]
678+
679+
if typeof(_f) <: ODEForwardSensitivityFunction
680+
_u0 = u0 === nothing ? state_values(prob, 1:(_f.numindvar)) :
681+
u0[1:(_f.numindvar)]
682+
else
683+
_u0 = u0 === nothing ? state_values(prob) : u0
684+
end
685+
680686
_tspan = tspan === nothing ? prob.tspan : tspan
681687
ODEForwardSensitivityProblem(_f, _u0,
682688
_tspan, _p; sensealg = prob.problem_type.sensealg,

test/forward.jl

Lines changed: 12 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -271,3 +271,15 @@ f = prob.f
271271
@assert f isa ODEForwardSensitivityFunction
272272
@test hasproperty(f, :observed)
273273
@test f.observed == SciMLBase.DEFAULT_OBSERVED
274+
275+
# `remake`: https://github.com/SciML/SciMLSensitivity.jl/issues/1137
276+
277+
function f(du, u, p, t)
278+
du[1] = dx = p[1] * u[1] - p[2] * u[1] * u[2]
279+
du[2] = dy = -p[3] * u[2] + u[1] * u[2]
280+
end
281+
282+
p = [1.5, 1.0, 3.0]
283+
ts = (0, 10)
284+
prob = ODEForwardSensitivityProblem(f, [1.0; 1.0], ts, p, sensealg=ForwardDiffSensitivity())
285+
sol = solve(prob, Tsit5())

0 commit comments

Comments
 (0)