/
DiffEqBaseReverseDiffExt.jl
173 lines (157 loc) · 6.47 KB
/
DiffEqBaseReverseDiffExt.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
module DiffEqBaseReverseDiffExt
if isdefined(Base, :get_extension)
using DiffEqBase
import DiffEqBase: value
import ReverseDiff
else
using ..DiffEqBase
import ..DiffEqBase: value
import ..ReverseDiff
end
DiffEqBase.value(x::Type{ReverseDiff.TrackedReal{V, D, O}}) where {V, D, O} = V
function DiffEqBase.value(x::Type{
ReverseDiff.TrackedArray{V, D, N, VA, DA},
}) where {V, D,
N, VA,
DA}
Array{V, N}
end
DiffEqBase.value(x::ReverseDiff.TrackedReal) = x.value
DiffEqBase.value(x::ReverseDiff.TrackedArray) = x.value
# Force TrackedArray from TrackedReal when reshaping W\b
DiffEqBase._reshape(v::AbstractVector{<:ReverseDiff.TrackedReal}, siz) = reduce(vcat, v)
DiffEqBase.promote_u0(u0::ReverseDiff.TrackedArray, p::ReverseDiff.TrackedArray, t0) = u0
function DiffEqBase.promote_u0(u0::AbstractArray{<:ReverseDiff.TrackedReal},
p::ReverseDiff.TrackedArray, t0)
u0
end
function DiffEqBase.promote_u0(u0::ReverseDiff.TrackedArray,
p::AbstractArray{<:ReverseDiff.TrackedReal}, t0)
u0
end
function DiffEqBase.promote_u0(u0::AbstractArray{<:ReverseDiff.TrackedReal},
p::AbstractArray{<:ReverseDiff.TrackedReal}, t0)
u0
end
DiffEqBase.promote_u0(u0, p::ReverseDiff.TrackedArray, t0) = ReverseDiff.track(u0)
DiffEqBase.promote_u0(u0, p::AbstractArray{<:ReverseDiff.TrackedReal}, t0) = eltype(p).(u0)
# Support adaptive with non-tracked time
@inline function DiffEqBase.ODE_DEFAULT_NORM(u::ReverseDiff.TrackedArray, t)
sqrt(sum(abs2, DiffEqBase.value(u)) / length(u))
end
@inline function DiffEqBase.ODE_DEFAULT_NORM(
u::AbstractArray{<:ReverseDiff.TrackedReal, N},
t) where {N}
sqrt(sum(x -> DiffEqBase.ODE_DEFAULT_NORM(x[1], x[2]),
zip((DiffEqBase.value(x) for x in u), Iterators.repeated(t))) / length(u))
end
@inline function DiffEqBase.ODE_DEFAULT_NORM(u::Array{<:ReverseDiff.TrackedReal, N},
t) where {N}
sqrt(sum(x -> DiffEqBase.ODE_DEFAULT_NORM(x[1], x[2]),
zip((DiffEqBase.value(x) for x in u), Iterators.repeated(t))) / length(u))
end
@inline function DiffEqBase.ODE_DEFAULT_NORM(u::ReverseDiff.TrackedReal, t)
abs(DiffEqBase.value(u))
end
# Support TrackedReal time, don't drop tracking on the adaptivity there
@inline function DiffEqBase.ODE_DEFAULT_NORM(u::ReverseDiff.TrackedArray,
t::ReverseDiff.TrackedReal)
sqrt(sum(abs2, u) / length(u))
end
@inline function DiffEqBase.ODE_DEFAULT_NORM(
u::AbstractArray{<:ReverseDiff.TrackedReal, N},
t::ReverseDiff.TrackedReal) where {N}
sqrt(sum(x -> DiffEqBase.ODE_DEFAULT_NORM(x[1], x[2]), zip(u, Iterators.repeated(t))) /
length(u))
end
@inline function DiffEqBase.ODE_DEFAULT_NORM(u::Array{<:ReverseDiff.TrackedReal, N},
t::ReverseDiff.TrackedReal) where {N}
sqrt(sum(x -> DiffEqBase.ODE_DEFAULT_NORM(x[1], x[2]), zip(u, Iterators.repeated(t))) /
length(u))
end
@inline function DiffEqBase.ODE_DEFAULT_NORM(u::ReverseDiff.TrackedReal,
t::ReverseDiff.TrackedReal)
abs(u)
end
# `ReverseDiff.TrackedArray`
function DiffEqBase.solve_up(prob::DiffEqBase.AbstractDEProblem,
sensealg::Union{
SciMLBase.AbstractOverloadingSensitivityAlgorithm,
Nothing}, u0::ReverseDiff.TrackedArray,
p::ReverseDiff.TrackedArray, args...; kwargs...)
ReverseDiff.track(DiffEqBase.solve_up, prob, sensealg, u0, p, args...; kwargs...)
end
function DiffEqBase.solve_up(prob::DiffEqBase.AbstractDEProblem,
sensealg::Union{
SciMLBase.AbstractOverloadingSensitivityAlgorithm,
Nothing}, u0, p::ReverseDiff.TrackedArray,
args...; kwargs...)
ReverseDiff.track(DiffEqBase.solve_up, prob, sensealg, u0, p, args...; kwargs...)
end
function DiffEqBase.solve_up(prob::DiffEqBase.AbstractDEProblem,
sensealg::Union{
SciMLBase.AbstractOverloadingSensitivityAlgorithm,
Nothing}, u0::ReverseDiff.TrackedArray, p,
args...; kwargs...)
ReverseDiff.track(DiffEqBase.solve_up, prob, sensealg, u0, p, args...; kwargs...)
end
# `AbstractArray{<:ReverseDiff.TrackedReal}`
function DiffEqBase.solve_up(prob::DiffEqBase.AbstractDEProblem,
sensealg::Union{
SciMLBase.AbstractOverloadingSensitivityAlgorithm,
Nothing},
u0::AbstractArray{<:ReverseDiff.TrackedReal},
p::AbstractArray{<:ReverseDiff.TrackedReal}, args...;
kwargs...)
DiffEqBase.solve_up(prob, sensealg, reduce(vcat, u0), reduce(vcat, p), args...;
kwargs...)
end
function DiffEqBase.solve_up(prob::DiffEqBase.AbstractDEProblem,
sensealg::Union{
SciMLBase.AbstractOverloadingSensitivityAlgorithm,
Nothing}, u0,
p::AbstractArray{<:ReverseDiff.TrackedReal},
args...; kwargs...)
DiffEqBase.solve_up(prob, sensealg, u0, reduce(vcat, p), args...; kwargs...)
end
function DiffEqBase.solve_up(prob::DiffEqBase.AbstractDEProblem,
sensealg::Union{
SciMLBase.AbstractOverloadingSensitivityAlgorithm,
Nothing}, u0::ReverseDiff.TrackedArray,
p::AbstractArray{<:ReverseDiff.TrackedReal},
args...; kwargs...)
DiffEqBase.solve_up(prob, sensealg, u0, reduce(vcat, p), args...; kwargs...)
end
function DiffEqBase.solve_up(prob::DiffEqBase.DEProblem,
sensealg::Union{
SciMLBase.AbstractOverloadingSensitivityAlgorithm,
Nothing},
u0::AbstractArray{<:ReverseDiff.TrackedReal}, p,
args...; kwargs...)
DiffEqBase.solve_up(prob, sensealg, reduce(vcat, u0), p, args...; kwargs...)
end
function DiffEqBase.solve_up(prob::DiffEqBase.DEProblem,
sensealg::Union{
SciMLBase.AbstractOverloadingSensitivityAlgorithm,
Nothing},
u0::AbstractArray{<:ReverseDiff.TrackedReal}, p::ReverseDiff.TrackedArray,
args...; kwargs...)
DiffEqBase.solve_up(prob, sensealg, reduce(vcat, u0), p, args...; kwargs...)
end
# Required becase ReverseDiff.@grad function DiffEqBase.solve_up is not supported!
import DiffEqBase: solve_up
ReverseDiff.@grad function solve_up(prob, sensealg, u0, p, args...; kwargs...)
out = DiffEqBase._solve_adjoint(prob, sensealg, ReverseDiff.value(u0),
ReverseDiff.value(p),
SciMLBase.ReverseDiffOriginator(), args...; kwargs...)
function actual_adjoint(_args...)
original_adjoint = out[2](_args...)
if isempty(args) # alg is missing
tuple(original_adjoint[1:4]..., original_adjoint[6:end]...)
else
original_adjoint
end
end
Array(out[1]), actual_adjoint
end
end