-
-
Notifications
You must be signed in to change notification settings - Fork 69
/
steady_state.jl
565 lines (504 loc) · 25.8 KB
/
steady_state.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
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
using Test, LinearAlgebra
using SciMLSensitivity, SteadyStateDiffEq, DiffEqBase, NLsolve
using OrdinaryDiffEq
using NonlinearSolve, SciMLNLSolve
using ForwardDiff, Calculus
using Zygote
using Random
Random.seed!(12345)
@testset "Adjoint sensitivities of steady state solver" begin
function f!(du, u, p, t)
du[1] = p[1] + p[2] * u[1]
du[2] = p[3] * u[1] + p[4] * u[2]
end
function jac!(J, u, p, t) #df/dx
J[1, 1] = p[2]
J[2, 1] = p[3]
J[1, 2] = 0
J[2, 2] = p[4]
nothing
end
function paramjac!(fp, u, p, t) #df/dp
fp[1, 1] = 1
fp[2, 1] = 0
fp[1, 2] = u[1]
fp[2, 2] = 0
fp[1, 3] = 0
fp[2, 3] = u[1]
fp[1, 4] = 0
fp[2, 4] = u[2]
nothing
end
function dgdu!(out, u, p, t, i)
(out .= -2.0 .+ u)
end
function dgdp!(out, u, p, t, i)
(out .= p)
end
function g(u, p, t)
sum((2.0 .- u) .^ 2) / 2 + sum(p .^ 2) / 2
end
u0 = zeros(2)
p = [2.0, -2.0, 1.0, -4.0]
prob = SteadyStateProblem(f!, u0, p)
abstol = 1e-10
@testset "for p" begin
println("Calculate adjoint sensitivities from Jacobians")
sol_analytical = [-p[1] / p[2], p[1] * p[3] / (p[2] * p[4])]
J = zeros(2, 2)
fp = zeros(2, 4)
gp = zeros(4)
gx = zeros(1, 2)
delg_delp = copy(p)
jac!(J, sol_analytical, p, nothing)
dgdu!(vec(gx), sol_analytical, p, nothing, nothing)
paramjac!(fp, sol_analytical, p, nothing)
lambda = J' \ gx'
res_analytical = delg_delp' - lambda' * fp # = -gx*inv(J)*fp
@info "Expected result" sol_analytical, res_analytical,
delg_delp' - gx * inv(J) * fp
@info "Calculate adjoint sensitivities from autodiff & numerical diff"
function G(p)
tmp_prob = remake(prob, u0 = convert.(eltype(p), prob.u0), p = p)
sol = solve(tmp_prob, DynamicSS(Rodas5()))
A = convert(Array, sol)
g(A, p, nothing)
end
res1 = ForwardDiff.gradient(G, p)
res2 = Calculus.gradient(G, p)
#@info res1, res2, res_analytical
@test res1≈res_analytical' rtol=1e-7
@test res2≈res_analytical' rtol=1e-7
@test res1≈res2 rtol=1e-7
@info "Adjoint sensitivities"
# with jac, param_jac
f1 = ODEFunction(f!; jac = jac!, paramjac = paramjac!)
prob1 = SteadyStateProblem(f1, u0, p)
sol1 = solve(prob1, DynamicSS(Rodas5(), reltol = 1e-14, abstol = 1e-14),
reltol = 1e-14, abstol = 1e-14)
res1a = adjoint_sensitivities(sol1, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), dgdu = dgdu!,
dgdp = dgdp!, g = g)
res1b = adjoint_sensitivities(sol1, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), g = g)
res1c = adjoint_sensitivities(sol1, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = false),
g = g)
res1d = adjoint_sensitivities(sol1, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = TrackerVJP()),
g = g)
res1e = adjoint_sensitivities(sol1, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = ReverseDiffVJP()),
g = g)
res1f = adjoint_sensitivities(sol1, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = ZygoteVJP()),
g = g)
res1g = adjoint_sensitivities(sol1, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = false,
autojacvec = false),
g = g)
res1h = adjoint_sensitivities(sol1, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = EnzymeVJP()),
g = g)
# with jac, without param_jac
f2 = ODEFunction(f!; jac = jac!)
prob2 = SteadyStateProblem(f2, u0, p)
sol2 = solve(prob2, DynamicSS(Rodas5(), reltol = 1e-14, abstol = 1e-14),
reltol = 1e-14, abstol = 1e-14)
res2a = adjoint_sensitivities(sol2, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), dgdu = dgdu!,
dgdp = dgdp!, g = g)
res2b = adjoint_sensitivities(sol2, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), g = g)
res2c = adjoint_sensitivities(sol2, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = false),
g = g)
res2d = adjoint_sensitivities(sol2, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = TrackerVJP()),
g = g)
res2e = adjoint_sensitivities(sol2, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = ReverseDiffVJP()),
g = g)
res2f = adjoint_sensitivities(sol2, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = ZygoteVJP()),
g = g)
res2g = adjoint_sensitivities(sol2, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = false,
autojacvec = false),
g = g)
res2h = adjoint_sensitivities(sol2, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = EnzymeVJP()),
g = g)
# without jac, without param_jac
f3 = ODEFunction(f!)
prob3 = SteadyStateProblem(f3, u0, p)
sol3 = solve(prob3, DynamicSS(Rodas5(), reltol = 1e-14, abstol = 1e-14),
reltol = 1e-14, abstol = 1e-14)
res3a = adjoint_sensitivities(sol3, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), dgdu = dgdu!,
dgdp = dgdp!, g = g)
res3b = adjoint_sensitivities(sol3, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), g = g)
res3c = adjoint_sensitivities(sol3, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = false),
g = g)
res3d = adjoint_sensitivities(sol3, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = TrackerVJP()),
g = g)
res3e = adjoint_sensitivities(sol3, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = ReverseDiffVJP()),
g = g)
res3f = adjoint_sensitivities(sol3, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = ZygoteVJP()),
g = g)
res3g = adjoint_sensitivities(sol3, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = false,
autojacvec = false),
g = g)
res3h = adjoint_sensitivities(sol3, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = EnzymeVJP()),
g = g)
@test norm(res_analytical' .- res1a) < 1e-7
@test norm(res_analytical' .- res1b) < 1e-7
@test norm(res_analytical' .- res1c) < 1e-7
@test norm(res_analytical' .- res1d) < 1e-7
@test norm(res_analytical' .- res1e) < 1e-7
@test norm(res_analytical' .- res1f) < 1e-7
@test norm(res_analytical' .- res1g) < 1e-7
@test norm(res_analytical' .- res1h) < 1e-7
@test norm(res_analytical' .- res2a) < 1e-7
@test norm(res_analytical' .- res2b) < 1e-7
@test norm(res_analytical' .- res2c) < 1e-7
@test norm(res_analytical' .- res2d) < 1e-7
@test norm(res_analytical' .- res2e) < 1e-7
@test norm(res_analytical' .- res2f) < 1e-7
@test norm(res_analytical' .- res2g) < 1e-7
@test norm(res_analytical' .- res2h) < 1e-7
@test norm(res_analytical' .- res3a) < 1e-7
@test norm(res_analytical' .- res3b) < 1e-7
@test norm(res_analytical' .- res3c) < 1e-7
@test norm(res_analytical' .- res3d) < 1e-7
@test norm(res_analytical' .- res3e) < 1e-7
@test norm(res_analytical' .- res3f) < 1e-7
@test norm(res_analytical' .- res3g) < 1e-7
@test norm(res_analytical' .- res3h) < 1e-7
@info "oop checks"
function foop(u, p, t)
dx = p[1] + p[2] * u[1]
dy = p[3] * u[1] + p[4] * u[2]
[dx, dy]
end
proboop = SteadyStateProblem(foop, u0, p)
soloop = solve(proboop, DynamicSS(Rodas5(), reltol = 1e-14, abstol = 1e-14),
reltol = 1e-14, abstol = 1e-14)
res4a = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), dgdu = dgdu!,
dgdp = dgdp!, g = g)
res4b = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), g = g)
res4c = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = false),
g = g)
res4d = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = TrackerVJP()),
g = g)
res4e = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = ReverseDiffVJP()),
g = g)
res4f = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autojacvec = ZygoteVJP()),
g = g)
res4g = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = false,
autojacvec = false),
g = g)
res4h = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(autodiff = true,
autojacvec = false),
g = g)
@test norm(res_analytical' .- res4a) < 1e-7
@test norm(res_analytical' .- res4b) < 1e-7
@test norm(res_analytical' .- res4c) < 1e-7
@test norm(res_analytical' .- res4d) < 1e-7
@test norm(res_analytical' .- res4e) < 1e-7
@test norm(res_analytical' .- res4f) < 1e-7
@test norm(res_analytical' .- res4g) < 1e-7
@test norm(res_analytical' .- res4h) < 1e-7
end
@testset "for u0: (should be zero, steady state does not depend on initial condition)" begin
res5 = ForwardDiff.gradient(prob.u0) do u0
tmp_prob = remake(prob, u0 = u0)
sol = solve(tmp_prob, DynamicSS(Rodas5()))
A = convert(Array, sol)
g(A, p, nothing)
end
@test abs(dot(res5, res5)) < 1e-7
end
end
@testset "concrete_solve derivatives steady state solver" begin
function g1(u, p, t)
sum(u)
end
function g2(u, p, t)
sum((2.0 .- u) .^ 2) / 2
end
u0 = zeros(2)
p = [2.0, -2.0, 1.0, -4.0]
@testset "iip" begin
function f!(du, u, p, t)
du[1] = p[1] + p[2] * u[1]
du[2] = p[3] * u[1] + p[4] * u[2]
end
prob = SteadyStateProblem(f!, u0, p)
sol = solve(prob, DynamicSS(Rodas5()))
res1 = adjoint_sensitivities(sol, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), g = g1)
res2 = adjoint_sensitivities(sol, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), g = g2)
dp1 = Zygote.gradient(p -> sum(solve(prob, DynamicSS(Rodas5()), u0 = u0, p = p,
sensealg = SteadyStateAdjoint())), p)
dp2 = Zygote.gradient(p -> sum((2.0 .-
solve(prob, DynamicSS(Rodas5()), u0 = u0, p = p,
sensealg = SteadyStateAdjoint())) .^ 2) / 2.0,
p)
dp1d = Zygote.gradient(p -> sum(solve(prob, DynamicSS(Rodas5()), u0 = u0, p = p)),
p)
dp2d = Zygote.gradient(p -> sum((2.0 .-
solve(prob, DynamicSS(Rodas5()), u0 = u0, p = p)) .^
2) / 2.0, p)
@test res1≈dp1[1] rtol=1e-12
@test res2≈dp2[1] rtol=1e-12
@test res1≈dp1d[1] rtol=1e-12
@test res2≈dp2d[1] rtol=1e-12
res1 = Zygote.gradient(p -> sum(Array(solve(prob, DynamicSS(Rodas5()), u0 = u0,
p = p, sensealg = SteadyStateAdjoint()))[1]),
p)
dp1 = Zygote.gradient(p -> sum(solve(prob, DynamicSS(Rodas5()), u0 = u0, p = p,
save_idxs = 1:1,
sensealg = SteadyStateAdjoint())), p)
dp2 = Zygote.gradient(p -> solve(prob, DynamicSS(Rodas5()), u0 = u0, p = p,
save_idxs = 1, sensealg = SteadyStateAdjoint())[1],
p)
dp1d = Zygote.gradient(p -> sum(solve(prob, DynamicSS(Rodas5()), u0 = u0, p = p,
save_idxs = 1:1)), p)
dp2d = Zygote.gradient(p -> solve(prob, DynamicSS(Rodas5()), u0 = u0, p = p,
save_idxs = 1)[1], p)
@test res1[1]≈dp1[1] rtol=1e-10
@test res1[1]≈dp2[1] rtol=1e-10
@test res1[1]≈dp1d[1] rtol=1e-10
@test res1[1]≈dp2d[1] rtol=1e-10
end
@testset "oop" begin
function f(u, p, t)
dx = p[1] + p[2] * u[1]
dy = p[3] * u[1] + p[4] * u[2]
[dx, dy]
end
proboop = SteadyStateProblem(f, u0, p)
soloop = solve(proboop, DynamicSS(Rodas5()))
res1oop = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), g = g1)
res2oop = adjoint_sensitivities(soloop, DynamicSS(Rodas5()),
sensealg = SteadyStateAdjoint(), g = g2)
dp1oop = Zygote.gradient(p -> sum(solve(proboop, DynamicSS(Rodas5()), u0 = u0,
p = p, sensealg = SteadyStateAdjoint())), p)
dp2oop = Zygote.gradient(p -> sum((2.0 .-
solve(proboop, DynamicSS(Rodas5()), u0 = u0,
p = p, sensealg = SteadyStateAdjoint())) .^
2) / 2.0, p)
dp1oopd = Zygote.gradient(p -> sum(solve(proboop, DynamicSS(Rodas5()), u0 = u0,
p = p)), p)
dp2oopd = Zygote.gradient(p -> sum((2.0 .-
solve(proboop, DynamicSS(Rodas5()), u0 = u0,
p = p)) .^ 2) / 2.0, p)
@test res1oop≈dp1oop[1] rtol=1e-12
@test res2oop≈dp2oop[1] rtol=1e-12
@test res1oop≈dp1oopd[1] rtol=1e-8
@test res2oop≈dp2oopd[1] rtol=1e-8
res1oop = Zygote.gradient(p -> sum(Array(solve(proboop, DynamicSS(Rodas5()),
u0 = u0, p = p,
sensealg = SteadyStateAdjoint()))[1]),
p)
dp1oop = Zygote.gradient(p -> sum(solve(proboop, DynamicSS(Rodas5()), u0 = u0,
p = p, save_idxs = 1:1,
sensealg = SteadyStateAdjoint())), p)
dp2oop = Zygote.gradient(p -> solve(proboop, DynamicSS(Rodas5()), u0 = u0, p = p,
save_idxs = 1, sensealg = SteadyStateAdjoint())[1],
p)
dp1oopd = Zygote.gradient(p -> sum(solve(proboop, DynamicSS(Rodas5()), u0 = u0,
p = p, save_idxs = 1:1)), p)
dp2oopd = Zygote.gradient(p -> solve(proboop, DynamicSS(Rodas5()), u0 = u0, p = p,
save_idxs = 1)[1], p)
@test res1oop[1]≈dp1oop[1] rtol=1e-10
@test res1oop[1]≈dp2oop[1] rtol=1e-10
end
end
@testset "NonlinearProblem" begin
u0 = [0.0]
p = [2.0, 1.0]
prob = NonlinearProblem((du, u, p) -> du[1] = u[1] - p[1] + p[2], u0, p)
prob2 = NonlinearProblem{false}((u, p) -> u .- p[1] .+ p[2], u0, p)
solve1 = solve(remake(prob, p = p), NewtonRaphson())
solve2 = solve(prob2, NewtonRaphson())
@test solve1.u == solve2.u
prob3 = SteadyStateProblem((u, p, t) -> -u .+ p[1] .- p[2], [0.0], p)
solve3 = solve(prob3, DynamicSS(Rodas5()))
@test solve1.u≈solve3.u rtol=1e-6
prob4 = SteadyStateProblem((du, u, p, t) -> du[1] = -u[1] + p[1] - p[2], [0.0], p)
solve4 = solve(prob4, DynamicSS(Rodas5()))
@test solve3.u≈solve4.u rtol=1e-10
function test_loss(p, prob; alg = NewtonRaphson())
_prob = remake(prob, p = p)
sol = sum(solve(_prob, alg,
sensealg = SteadyStateAdjoint(autojacvec = ReverseDiffVJP())))
return sol
end
test_loss(p, prob)
test_loss(p, prob2)
test_loss(p, prob3, alg = DynamicSS(Rodas5()))
test_loss(p, prob4, alg = DynamicSS(Rodas5()))
test_loss(p, prob2, alg = SimpleNewtonRaphson())
dp1 = Zygote.gradient(p -> test_loss(p, prob), p)[1]
dp2 = Zygote.gradient(p -> test_loss(p, prob2), p)[1]
dp3 = Zygote.gradient(p -> test_loss(p, prob3, alg = DynamicSS(Rodas5())), p)[1]
dp4 = Zygote.gradient(p -> test_loss(p, prob4, alg = DynamicSS(Rodas5())), p)[1]
dp5 = Zygote.gradient(p -> test_loss(p, prob2, alg = SimpleNewtonRaphson()), p)[1]
dp6 = Zygote.gradient(p -> test_loss(p, prob2, alg = Klement()), p)[1]
dp7 = Zygote.gradient(p -> test_loss(p, prob2, alg = SimpleTrustRegion()), p)[1]
dp8 = Zygote.gradient(p -> test_loss(p, prob2, alg = NLSolveJL()), p)[1]
@test dp1≈dp2 rtol=1e-10
@test dp1≈dp3 rtol=1e-10
@test dp1≈dp4 rtol=1e-10
@test dp1≈dp5 rtol=1e-10
@test dp1≈dp6 rtol=1e-10
@test dp1≈dp7 rtol=1e-10
@test dp1≈dp8 rtol=1e-10
end
@testset "Continuous sensitivity tools" begin
function f!(du, u, p, t)
du[1] = p[1] + p[2] * u[1]
du[2] = p[3] * u[1] + p[4] * u[2]
end
function g(u, p, t)
sum((2.0 .- u) .^ 2) / 2 + sum(p .^ 2) / 2
end
u0 = zeros(2)
p = [2.0, -2.0, 1.0, -4.0]
prob = ODEProblem(f!, u0, (0, Inf), p)
tol = 1e-10
# steady state callback
function condition(u, t, integrator)
testval = first(get_tmp_cache(integrator))
DiffEqBase.get_du!(testval, integrator)
all(testval .< tol)
end
affect!(integrator) = terminate!(integrator)
cb_t1 = DiscreteCallback(condition, affect!, save_positions = (false, true))
cb_t2 = DiscreteCallback(condition, affect!, save_positions = (true, true))
for cb_t in (cb_t1, cb_t2)
sol = solve(prob, Tsit5(), reltol = tol, abstol = tol, callback = cb_t,
save_start = false, save_everystep = false)
# derivative with respect to u0 and p0
function loss(u0, p; sensealg = nothing, save_start = false, save_everystep = false)
_prob = remake(prob, u0 = u0, p = p)
# saving arguments can have a huge influence here
sol = solve(_prob, Tsit5(), reltol = tol, abstol = tol, sensealg = sensealg,
callback = cb_t,
save_start = save_start, save_everystep = save_everystep)
res = sol.u[end]
g(res, p, nothing)
end
du0 = ForwardDiff.gradient((u0) -> loss(u0, p), u0)
dp = ForwardDiff.gradient((p) -> loss(u0, p), p)
# save_start = false, save_everystep=false
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p,
sensealg = ForwardDiffSensitivity()),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p, sensealg = BacksolveAdjoint()),
u0,
p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p,
sensealg = InterpolatingAdjoint()),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
# save_start = true, save_everystep=false
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p,
sensealg = ForwardDiffSensitivity(),
save_start = true), u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p, sensealg = BacksolveAdjoint(),
save_start = true), u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p,
sensealg = InterpolatingAdjoint(),
save_start = true), u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
# save_start = true, save_everystep=true
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p,
sensealg = ForwardDiffSensitivity(),
save_start = true,
save_everystep = true),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p, sensealg = BacksolveAdjoint(),
save_start = true,
save_everystep = true),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p,
sensealg = InterpolatingAdjoint(),
save_start = true,
save_everystep = true),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
# QuadratureAdjoint makes sense only in this case, otherwise Zdp fails
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss(u0, p, sensealg = QuadratureAdjoint(),
save_start = true,
save_everystep = true),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
function loss2(u0, p; sensealg = nothing, saveat = 1.0)
# remake tspan so saveat::Number makes sense
_prob = remake(prob, tspan = (0.0, 100.0), u0 = u0, p = p)
# saving arguments can have a huge influence here
sol = solve(_prob, Tsit5(), reltol = tol, abstol = tol, sensealg = sensealg,
callback = cb_t, saveat = saveat)
res = sol.u[end]
g(res, p, nothing)
end
du0 = ForwardDiff.gradient((u0) -> loss2(u0, p), u0)
dp = ForwardDiff.gradient((p) -> loss2(u0, p), p)
# saveat::Number
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss2(u0, p,
sensealg = ForwardDiffSensitivity()),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss2(u0, p, sensealg = BacksolveAdjoint()),
u0,
p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss2(u0, p,
sensealg = InterpolatingAdjoint()),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
Zdu0, Zdp = Zygote.gradient((u0, p) -> loss2(u0, p, sensealg = QuadratureAdjoint()),
u0, p)
@test du0≈Zdu0 atol=1e-4
@test dp≈Zdp atol=1e-4
end
end