-
-
Notifications
You must be signed in to change notification settings - Fork 89
/
SciMLBaseZygoteExt.jl
308 lines (278 loc) · 11 KB
/
SciMLBaseZygoteExt.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
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
module SciMLBaseZygoteExt
using Zygote
using Zygote: @adjoint, pullback
import Zygote: literal_getproperty
using SciMLBase
using SciMLBase: ODESolution, remake,
getobserved, build_solution, EnsembleSolution,
NonlinearSolution, AbstractTimeseriesSolution
using SymbolicIndexingInterface: symbolic_type, NotSymbolic, variable_index, is_observed,
observed, parameter_values, state_values, current_time
using RecursiveArrayTools
import SciMLStructures
# This method resolves the ambiguity with the pullback defined in
# RecursiveArrayToolsZygoteExt
# https://github.com/SciML/RecursiveArrayTools.jl/blob/d06ecb856f43bc5e37cbaf50e5f63c578bf3f1bd/ext/RecursiveArrayToolsZygoteExt.jl#L67
@adjoint function Base.getindex(VA::ODESolution, i::Int, j::Int)
function ODESolution_getindex_pullback(Δ)
du = [m == j ? [i == k ? Δ : zero(VA.u[1][1]) for k in 1:length(VA.u[1])] :
zero(VA.u[1]) for m in 1:length(VA.u)]
dp = zero(VA.prob.p)
dprob = remake(VA.prob, p = dp)
du, dprob
T = eltype(eltype(VA.u))
if dprob.u0 === nothing
N = 2
elseif dprob isa SciMLBase.BVProblem && !hasmethod(size, Tuple{typeof(dprob.u0)})
__u0 = hasmethod(dprob.u0, Tuple{typeof(dprob.p), typeof(first(dprob.tspan))}) ?
dprob.u0(dprob.p, first(dprob.tspan)) : dprob.u0(first(dprob.tspan))
N = length((size(__u0)..., length(du)))
else
N = length((size(dprob.u0)..., length(du)))
end
Δ′ = ODESolution{T, N}(du, nothing, nothing,
VA.t, VA.k, dprob, VA.alg, VA.interp, VA.dense, 0, VA.stats,
VA.alg_choice, VA.retcode)
(Δ′, nothing, nothing)
end
VA[i, j], ODESolution_getindex_pullback
end
@adjoint function Base.getindex(VA::ODESolution, sym, j::Int)
function ODESolution_getindex_pullback(Δ)
i = symbolic_type(sym) != NotSymbolic() ? variable_index(VA, sym) : sym
du, dprob = if i === nothing
getter = getobserved(VA)
grz = pullback(getter, sym, VA.u[j], VA.prob.p, VA.t[j])[2](Δ)
du = [k == j ? grz[2] : zero(VA.u[1]) for k in 1:length(VA.u)]
dp = grz[3] # pullback for p
dprob = remake(VA.prob, p = dp)
du, dprob
else
du = [m == j ? [i == k ? Δ : zero(VA.u[1][1]) for k in 1:length(VA.u[1])] :
zero(VA.u[1]) for m in 1:length(VA.u)]
dp = zero(VA.prob.p)
dprob = remake(VA.prob, p = dp)
du, dprob
end
T = eltype(eltype(VA.u))
if dprob.u0 === nothing
N = 2
elseif dprob isa SciMLBase.BVProblem && !hasmethod(size, Tuple{typeof(dprob.u0)})
__u0 = hasmethod(dprob.u0, Tuple{typeof(dprob.p), typeof(first(dprob.tspan))}) ?
dprob.u0(dprob.p, first(dprob.tspan)) : dprob.u0(first(dprob.tspan))
N = length((size(__u0)..., length(du)))
else
N = length((size(dprob.u0)..., length(du)))
end
Δ′ = ODESolution{T, N}(du, nothing, nothing,
VA.t, VA.k, dprob, VA.alg, VA.interp, VA.dense, 0, VA.stats,
VA.alg_choice, VA.retcode)
(Δ′, nothing, nothing)
end
VA[sym, j], ODESolution_getindex_pullback
end
@adjoint function EnsembleSolution(sim, time, converged, stats)
out = EnsembleSolution(sim, time, converged, stats)
function EnsembleSolution_adjoint(p̄::AbstractArray{T, N}) where {T, N}
arrarr = [[p̄[ntuple(x -> Colon(), Val(N - 2))..., j, i]
for j in 1:size(p̄)[end - 1]] for i in 1:size(p̄)[end]]
(EnsembleSolution(arrarr, 0.0, true, stats), nothing, nothing, nothing)
end
function EnsembleSolution_adjoint(p̄::AbstractArray{<:AbstractArray, 1})
(EnsembleSolution(p̄, 0.0, true, stats), nothing, nothing, nothing)
end
function EnsembleSolution_adjoint(p̄::RecursiveArrayTools.AbstractVectorOfArray)
(EnsembleSolution(p̄, 0.0, true, stats), nothing, nothing, nothing)
end
function EnsembleSolution_adjoint(p̄::EnsembleSolution)
(p̄, nothing, nothing, nothing)
end
out, EnsembleSolution_adjoint
end
@adjoint function Base.getindex(VA::ODESolution, i::Int)
function ODESolution_getindex_pullback(Δ)
Δ′ = [(i == j ? Δ : Zygote.FillArrays.Fill(zero(eltype(x)), size(x)))
for (x, j) in zip(VA.u, 1:length(VA))]
(Δ′, nothing)
end
VA[:, i], ODESolution_getindex_pullback
end
@adjoint function Zygote.literal_getproperty(sim::EnsembleSolution,
::Val{:u})
sim.u, p̄ -> (EnsembleSolution(p̄, 0.0, true, sim.stats),)
end
@adjoint function Base.getindex(VA::ODESolution, sym)
function ODESolution_getindex_pullback(Δ)
i = symbolic_type(sym) != NotSymbolic() ? variable_index(VA, sym) : sym
if is_observed(VA, sym)
f = observed(VA, sym)
p = parameter_values(VA)
tunables, _, _ = SciMLStructures.canonicalize(SciMLStructures.Tunable(), p)
u = state_values(VA)
t = current_time(VA)
y, back = Zygote.pullback(u, tunables) do u, tunables
f.(u, Ref(tunables), t)
end
gs = back(Δ)
(gs[1], nothing)
elseif i === nothing
throw(error("Zygote AD of purely-symbolic slicing for observed quantities is not yet supported. Work around this by using `A[sym,i]` to access each element sequentially in the function being differentiated."))
else
Δ′ = [[i == k ? Δ[j] : zero(x[1]) for k in 1:length(x)]
for (x, j) in zip(VA.u, 1:length(VA))]
(Δ′, nothing)
end
end
VA[sym], ODESolution_getindex_pullback
end
function obs_grads(VA, sym, obs_idx, Δ)
y, back = Zygote.pullback(VA) do sol
getindex.(Ref(sol), sym[obs_idx])
end
Δreduced = reduce(hcat, Δ)
Δobs = eachrow(Δreduced[obs_idx, :])
back(Δobs)
end
function obs_grads(VA, sym, ::Nothing, Δ)
Zygote.nt_nothing(VA)
end
function not_obs_grads(VA::ODESolution{T}, sym, not_obss_idx, i, Δ) where {T}
Δ′ = map(enumerate(VA.u)) do (t_idx, us)
map(enumerate(us)) do (u_idx, u)
if u_idx in i
idx = findfirst(isequal(u_idx), i)
Δ[t_idx][idx]
else
zero(T)
end
end
end
Δ′
end
@adjoint function Base.getindex(
VA::ODESolution{T}, sym::Union{Tuple, AbstractVector}) where {T}
function ODESolution_getindex_pullback(Δ)
sym = sym isa Tuple ? collect(sym) : sym
i = map(x -> symbolic_type(x) != NotSymbolic() ? variable_index(VA, x) : x, sym)
obs_idx = findall(s -> is_observed(VA, s), sym)
not_obs_idx = setdiff(1:length(sym), obs_idx)
gs_obs = obs_grads(VA, sym, isempty(obs_idx) ? nothing : obs_idx, Δ)
gs_not_obs = not_obs_grads(VA, sym, not_obs_idx, i, Δ)
a = Zygote.accum(gs_obs[1], gs_not_obs)
(a, nothing)
end
VA[sym], ODESolution_getindex_pullback
end
@adjoint function ODESolution{T1, T2, T3, T4, T5, T6, T7, T8, T9, T10, T11, T12
}(u,
args...) where {T1, T2, T3, T4, T5, T6, T7, T8,
T9, T10, T11, T12}
function ODESolutionAdjoint(ȳ)
(ȳ, ntuple(_ -> nothing, length(args))...)
end
ODESolution{T1, T2, T3, T4, T5, T6, T7, T8, T9, T10, T11, T12}(u, args...),
ODESolutionAdjoint
end
@adjoint function SDEProblem{uType, tType, isinplace, P, NP, F, G, K, ND}(u,
args...) where
{uType, tType, isinplace, P, NP, F, G, K, ND}
function SDEProblemAdjoint(ȳ)
(ȳ, ntuple(_ -> nothing, length(args))...)
end
SDEProblem{uType, tType, isinplace, P, NP, F, G, K, ND}(u, args...), SDEProblemAdjoint
end
@adjoint function NonlinearSolution{T, N, uType, R, P, A, O, uType2}(u,
args...) where {
T,
N,
uType,
R,
P,
A,
O,
uType2
}
function NonlinearSolutionAdjoint(ȳ)
(ȳ, ntuple(_ -> nothing, length(args))...)
end
NonlinearSolution{T, N, uType, R, P, A, O, uType2}(u, args...), NonlinearSolutionAdjoint
end
@adjoint function literal_getproperty(sol::AbstractTimeseriesSolution,
::Val{:u})
function solu_adjoint(Δ)
zerou = zero(sol.prob.u0)
_Δ = @. ifelse(Δ === nothing, (zerou,), Δ)
(build_solution(sol.prob, sol.alg, sol.t, _Δ),)
end
sol.u, solu_adjoint
end
@adjoint function literal_getproperty(sol::SciMLBase.AbstractNoTimeSolution,
::Val{:u})
function solu_adjoint(Δ)
zerou = zero(sol.prob.u0)
_Δ = @. ifelse(Δ === nothing, zerou, Δ)
(build_solution(sol.prob, sol.alg, _Δ, sol.resid),)
end
sol.u, solu_adjoint
end
@adjoint function literal_getproperty(sol::SciMLBase.LinearSolution, ::Val{:u})
function solu_adjoint(Δ)
zerou = zero(sol.u)
_Δ = @. ifelse(Δ === nothing, zerou, Δ)
(SciMLBase.build_linear_solution(sol.cache.alg, _Δ, sol.resid, sol.cache),)
end
sol.u, solu_adjoint
end
@adjoint function literal_getproperty(sol::SciMLBase.OptimizationSolution,
::Val{:u})
function solu_adjoint(Δ)
zerou = zero(sol.u)
_Δ = @. ifelse(Δ === nothing, zerou, Δ)
(build_solution(sol.cache, sol.alg, _Δ, sol.objective),)
end
sol.u, solu_adjoint
end
function ∇tmap(cx, f, args...)
ys_and_backs = SciMLBase.tmap((args...) -> Zygote._pullback(cx, f, args...), args...)
if isempty(ys_and_backs)
ys_and_backs, _ -> (NoTangent(), NoTangent())
else
ys, backs = Zygote.unzip(ys_and_backs)
function ∇tmap_internal(Δ)
Δf_and_args_zipped = SciMLBase.tmap((f, δ) -> f(δ), backs, Δ)
Δf_and_args = Zygote.unzip(Δf_and_args_zipped)
Δf = reduce(Zygote.accum, Δf_and_args[1])
(Δf, Δf_and_args[2:end]...)
end
ys, ∇tmap_internal
end
end
function ∇responsible_map(cx, f, args...)
ys_and_backs = SciMLBase.responsible_map((args...) -> Zygote._pullback(cx, f, args...),
args...)
if isempty(ys_and_backs)
ys_and_backs, _ -> (NoTangent(), NoTangent())
else
ys, backs = Zygote.unzip(ys_and_backs)
ys,
function ∇responsible_map_internal(Δ)
# Apply pullbacks in reverse order. Needed for correctness if `f` is stateful.
Δf_and_args_zipped = SciMLBase.responsible_map((f, δ) -> f(δ),
Zygote._tryreverse(SciMLBase.responsible_map,
backs, Δ)...)
Δf_and_args = Zygote.unzip(Zygote._tryreverse(SciMLBase.responsible_map,
Δf_and_args_zipped))
Δf = reduce(Zygote.accum, Δf_and_args[1])
(Δf, Δf_and_args[2:end]...)
end
end
end
@adjoint function SciMLBase.tmap(f, args::Union{AbstractArray, Tuple}...)
∇tmap(__context__, f, args...)
end
@adjoint function SciMLBase.responsible_map(f,
args::Union{AbstractArray, Tuple
}...)
∇responsible_map(__context__, f, args...)
end
end