/
1_basic.jl
242 lines (181 loc) · 5.93 KB
/
1_basic.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
# # Basic use cases
#=
We show how to differentiate through very common routines:
- an unconstrained optimization problem
- a nonlinear system of equations
- a fixed point iteration
=#
using ForwardDiff
using ImplicitDifferentiation
using LinearAlgebra
using NLsolve
using Optim
using Random
using Test #src
using Zygote
Random.seed!(63);
#=
In all three cases, we will use the square root as our forward mapping, but expressed in three different ways.
Here's our heroic test vector:
=#
x = rand(2);
#=
Since we already know the mathematical expression of the Jacobian, we will be able to compare it with our numerical results.
=#
J = Diagonal(0.5 ./ sqrt.(x))
# ## Unconstrained optimization
#=
First, we show how to differentiate through the solution of an unconstrained optimization problem:
```math
y(x) = \underset{y \in \mathbb{R}^m}{\mathrm{argmin}} ~ f(x, y)
```
The optimality conditions are given by gradient stationarity:
```math
\nabla_2 f(x, y) = 0
```
=#
#=
To make verification easy, we minimize the following objective:
```math
f(x, y) = \lVert y \odot y - x \rVert^2
```
In this case, the optimization problem boils down to the componentwise square root function, but we implement it using a black box solver from [Optim.jl](https://github.com/JuliaNLSolvers/Optim.jl).
Note the presence of a keyword argument.
=#
function forward_optim(x; method)
f(y) = sum(abs2, y .^ 2 .- x)
y0 = ones(eltype(x), size(x))
result = optimize(f, y0, method)
return Optim.minimizer(result)
end
#=
Even though they are defined as a gradient, it is better to provide optimality conditions explicitly: that way we avoid nesting autodiff calls. By default, the conditions should accept two arguments as input.
The forward mapping and the conditions should accept the same set of keyword arguments.
=#
function conditions_optim(x, y; method)
∇₂f = @. 4 * (y^2 - x) * y
return ∇₂f
end
#=
We now have all the ingredients to construct our implicit function.
=#
implicit_optim = ImplicitFunction(forward_optim, conditions_optim)
# And indeed, it behaves as it should when we call it:
implicit_optim(x; method=LBFGS()) .^ 2
@test implicit_optim(x; method=LBFGS()) .^ 2 ≈ x #src
# Forward mode autodiff
ForwardDiff.jacobian(_x -> implicit_optim(_x; method=LBFGS()), x)
@test ForwardDiff.jacobian(_x -> implicit_optim(_x; method=LBFGS()), x) ≈ J #src
#=
In this instance, we could use ForwardDiff.jl directly on the solver, but it returns the wrong result (not sure why).
=#
ForwardDiff.jacobian(_x -> forward_optim(x; method=LBFGS()), x)
# Reverse mode autodiff
Zygote.jacobian(_x -> implicit_optim(_x; method=LBFGS()), x)[1]
@test Zygote.jacobian(_x -> implicit_optim(_x; method=LBFGS()), x)[1] ≈ J #src
#=
In this instance, we cannot use Zygote.jl directly on the solver (due to unsupported `try/catch` statements).
=#
try
Zygote.jacobian(_x -> forward_optim(x; method=LBFGS()), x)[1]
catch e
e
end
# ## Nonlinear system
#=
Next, we show how to differentiate through the solution of a nonlinear system of equations:
```math
\text{find} \quad y(x) \quad \text{such that} \quad F(x, y(x)) = 0
```
The optimality conditions are pretty obvious:
```math
F(x, y) = 0
```
=#
#=
To make verification easy, we solve the following system:
```math
F(x, y) = y \odot y - x = 0
```
In this case, the optimization problem boils down to the componentwise square root function, but we implement it using a black box solver from [NLsolve.jl](https://github.com/JuliaNLSolvers/NLsolve.jl).
=#
function forward_nlsolve(x; method)
F!(storage, y) = (storage .= y .^ 2 .- x)
initial_y = similar(x)
initial_y .= 1
result = nlsolve(F!, initial_y; method)
return result.zero
end
#-
function conditions_nlsolve(x, y; method)
c = y .^ 2 .- x
return c
end
#-
implicit_nlsolve = ImplicitFunction(forward_nlsolve, conditions_nlsolve)
#-
implicit_nlsolve(x; method=:newton) .^ 2
@test implicit_nlsolve(x; method=:newton) .^ 2 ≈ x #src
# Forward mode autodiff
ForwardDiff.jacobian(_x -> implicit_nlsolve(_x; method=:newton), x)
@test ForwardDiff.jacobian(_x -> implicit_nlsolve(_x; method=:newton), x) ≈ J #src
#-
ForwardDiff.jacobian(_x -> forward_nlsolve(_x; method=:newton), x)
# Reverse mode autodiff
Zygote.jacobian(_x -> implicit_nlsolve(_x; method=:newton), x)[1]
@test Zygote.jacobian(_x -> implicit_nlsolve(_x; method=:newton), x)[1] ≈ J #src
#-
try
Zygote.jacobian(_x -> forward_nlsolve(_x; method=:newton), x)[1]
catch e
e
end
# ## Fixed point
#=
Finally, we show how to differentiate through the limit of a fixed point iteration:
```math
y \longmapsto T(x, y)
```
The optimality conditions are pretty obvious:
```math
y = T(x, y)
```
=#
#=
To make verification easy, we consider [Heron's method](https://en.wikipedia.org/wiki/Methods_of_computing_square_roots#Heron's_method):
```math
T(x, y) = \frac{1}{2} \left(y + \frac{x}{y}\right)
```
In this case, the fixed point algorithm boils down to the componentwise square root function, but we implement it manually.
=#
function forward_fixedpoint(x; iterations)
y = ones(eltype(x), size(x))
for _ in 1:iterations
y .= 0.5 .* (y .+ x ./ y)
end
return y
end
#-
function conditions_fixedpoint(x, y; iterations)
T = 0.5 .* (y .+ x ./ y)
return T .- y
end
#-
implicit_fixedpoint = ImplicitFunction(forward_fixedpoint, conditions_fixedpoint)
#-
implicit_fixedpoint(x; iterations=10) .^ 2
@test implicit_fixedpoint(x; iterations=10) .^ 2 ≈ x #src
# Forward mode autodiff
ForwardDiff.jacobian(_x -> implicit_fixedpoint(_x; iterations=10), x)
@test ForwardDiff.jacobian(_x -> implicit_fixedpoint(_x; iterations=10), x) ≈ J #src
#-
ForwardDiff.jacobian(_x -> forward_fixedpoint(_x; iterations=10), x)
# Reverse mode autodiff
Zygote.jacobian(_x -> implicit_fixedpoint(_x; iterations=10), x)[1]
@test Zygote.jacobian(_x -> implicit_fixedpoint(_x; iterations=10), x)[1] ≈ J #src
#-
try
Zygote.jacobian(_x -> forward_fixedpoint(_x; iterations=10), x)[1]
catch e
e
end