-
-
Notifications
You must be signed in to change notification settings - Fork 74
/
Copy pathcontinuous_vs_discrete.jl
193 lines (160 loc) · 5.97 KB
/
continuous_vs_discrete.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
using OrdinaryDiffEq, Zygote
using SciMLSensitivity, Test, ForwardDiff
abstol = 1e-12
reltol = 1e-12
savingtimes = 0.5
function test_continuous_wrt_discrete_callback()
# test the continuous callbacks wrt to the equivalent discrete callback
function f(du, u, p, t)
#Bouncing Ball
du[1] = u[2]
du[2] = -p[1]
end
# no saving in Callbacks; prescribed vafter and vbefore; loss on the endpoint
tstop = 3.1943828249997
vbefore = -31.30495168499705
vafter = 25.04396134799764
u0 = [50.0, 0.0]
tspan = (0.0, 5.0)
p = [9.8, 0.8]
prob = ODEProblem(f, u0, tspan, p)
function condition(u, t, integrator) # Event when event_f(u,t) == 0
t - tstop
end
function affect!(integrator)
integrator.u[2] += vafter - vbefore
end
cb = ContinuousCallback(condition, affect!, save_positions = (false, false))
condition2(u, t, integrator) = t == tstop
cb2 = DiscreteCallback(condition2, affect!, save_positions = (false, false))
du01, dp1 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb2, tstops = [tstop],
sensealg = BacksolveAdjoint(),
saveat = tspan[2], save_start = false)),
u0, p)
du02, dp2 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb,
sensealg = BacksolveAdjoint(),
saveat = tspan[2], save_start = false)),
u0, p)
dstuff = ForwardDiff.gradient(
(θ) -> sum(solve(prob, Tsit5(), u0 = θ[1:2], p = θ[3:4],
callback = cb, saveat = tspan[2],
save_start = false)),
[u0; p])
@info dstuff
@test du01 ≈ dstuff[1:2]
@test dp1 ≈ dstuff[3:4]
@test du02 ≈ dstuff[1:2]
@test dp2 ≈ dstuff[3:4]
# no saving in Callbacks; prescribed vafter and vbefore; loss on the endpoint by slicing
du01, dp1 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb2, tstops = [tstop],
sensealg = BacksolveAdjoint())[end]),
u0, p)
du02, dp2 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb,
sensealg = BacksolveAdjoint())[end]),
u0, p)
dstuff = ForwardDiff.gradient(
(θ) -> sum(solve(prob, Tsit5(), u0 = θ[1:2], p = θ[3:4],
callback = cb)[end]),
[u0; p])
@info dstuff
@test du01 ≈ dstuff[1:2]
@test dp1 ≈ dstuff[3:4]
@test du02 ≈ dstuff[1:2]
@test dp2 ≈ dstuff[3:4]
# with saving in Callbacks; prescribed vafter and vbefore; loss on the endpoint
cb = ContinuousCallback(condition, affect!, save_positions = (true, true))
cb2 = DiscreteCallback(condition2, affect!, save_positions = (true, true))
du01, dp1 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb2, tstops = [tstop],
sensealg = BacksolveAdjoint(),
saveat = tspan[2], save_start = false)),
u0, p)
du02, dp2 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb,
sensealg = BacksolveAdjoint(),
saveat = tspan[2], save_start = false)),
u0, p)
dstuff = ForwardDiff.gradient(
(θ) -> sum(solve(prob, Tsit5(), u0 = θ[1:2], p = θ[3:4],
callback = cb, saveat = tspan[2],
save_start = false)),
[u0; p])
@info dstuff
@test du01 ≈ dstuff[1:2]
@test dp1 ≈ dstuff[3:4]
@test du02 ≈ dstuff[1:2]
@test dp2 ≈ dstuff[3:4]
# with saving in Callbacks; prescribed vafter and vbefore; loss on the endpoint by slicing
du01, dp1 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb2, tstops = [tstop],
sensealg = BacksolveAdjoint())[end]),
u0, p)
du02, dp2 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb,
sensealg = BacksolveAdjoint())[end]),
u0, p)
dstuff = ForwardDiff.gradient(
(θ) -> sum(solve(prob, Tsit5(), u0 = θ[1:2], p = θ[3:4],
callback = cb)[end]),
[u0; p])
@info dstuff
@test du01 ≈ dstuff[1:2]
@test dp1 ≈ dstuff[3:4]
@test du02 ≈ dstuff[1:2]
@test dp2 ≈ dstuff[3:4]
# with saving in Callbacks; different affect function
function affect2!(integrator)
integrator.u[2] = -integrator.p[2] * integrator.u[2]
end
cb = ContinuousCallback(condition, affect2!, save_positions = (true, true))
cb2 = DiscreteCallback(condition2, affect2!, save_positions = (true, true))
du01, dp1 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb2, tstops = [tstop],
sensealg = BacksolveAdjoint(),
saveat = tspan[2], save_start = false)),
u0, p)
du02, dp2 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb,
sensealg = BacksolveAdjoint(),
saveat = tspan[2], save_start = false)),
u0, p)
du03, dp3 = Zygote.gradient(
(u0, p) -> sum(solve(prob, Tsit5(), u0 = u0, p = p,
callback = cb,
sensealg = GaussAdjoint(),
saveat = tspan[2], save_start = false)),
u0, p)
dstuff = ForwardDiff.gradient(
(θ) -> sum(solve(prob, Tsit5(), u0 = θ[1:2], p = θ[3:4],
callback = cb, saveat = tspan[2],
save_start = false)),
[u0; p])
@info dstuff
@test du01 ≈ dstuff[1:2]
@test dp1 ≈ dstuff[3:4]
@test du02 ≈ dstuff[1:2]
@test dp2 ≈ dstuff[3:4]
@test du03 ≈ dstuff[1:2]
@test dp3 ≈ dstuff[3:4]
@test du01 ≈ du02
@test dp1 ≈ dp2
@test du01 ≈ du03
@test dp1 ≈ dp3
end
@testset "Compare continuous with discrete callbacks" begin
test_continuous_wrt_discrete_callback()
end