-
-
Notifications
You must be signed in to change notification settings - Fork 203
/
utils.jl
412 lines (333 loc) · 12.6 KB
/
utils.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
get_status(nlsolver::AbstractNLSolver) = nlsolver.status
get_new_W_γdt_cutoff(nlsolver::AbstractNLSolver) = nlsolver.cache.new_W_γdt_cutoff
# handle FIRK
get_new_W_γdt_cutoff(alg::NewtonAlgorithm) = alg.new_W_γdt_cutoff
nlsolvefail(nlsolver::AbstractNLSolver) = nlsolvefail(get_status(nlsolver))
nlsolvefail(status::NLStatus) = Int8(status) <= 0
isnewton(::Any) = false
isnewton(nlsolver::AbstractNLSolver) = isnewton(nlsolver.cache)
isnewton(::AbstractNLSolverCache) = false
isnewton(::Union{NLNewtonCache,NLNewtonConstantCache}) = true
isJcurrent(nlsolver::AbstractNLSolver, integrator) = integrator.t == nlsolver.cache.J_t
isfirstcall(nlsolver::AbstractNLSolver) = nlsolver.cache.firstcall
isfirststage(nlsolver::AbstractNLSolver) = nlsolver.cache.firststage
setfirststage!(nlsolver::AbstractNLSolver, val::Bool) = setfirststage!(nlsolver.cache, val)
setfirststage!(nlcache::Union{NLNewtonCache,NLNewtonConstantCache}, val::Bool) = (nlcache.firststage = val)
setfirststage!(::Any, val::Bool) = nothing
markfirststage!(nlsolver::AbstractNLSolver) = setfirststage!(nlsolver, true)
set_new_W!(nlsolver::AbstractNLSolver, val::Bool)::Bool = set_new_W!(nlsolver.cache, val)
set_new_W!(nlcache::Union{NLNewtonCache,NLNewtonConstantCache}, val::Bool)::Bool =
nlcache.new_W = val
get_new_W!(nlsolver::AbstractNLSolver)::Bool = get_new_W!(nlsolver.cache)
get_new_W!(nlcache::Union{NLNewtonCache,NLNewtonConstantCache})::Bool = nlcache.new_W
get_new_W!(::AbstractNLSolverCache)::Bool = true
get_W(nlsolver::AbstractNLSolver) = get_W(nlsolver.cache)
get_W(nlcache::Union{NLNewtonCache,NLNewtonConstantCache}) = nlcache.W
set_W_γdt!(nlsolver::AbstractNLSolver, W_γdt) = set_W_γdt!(nlsolver.cache, W_γdt)
function set_W_γdt!(nlcache::Union{NLNewtonCache,NLNewtonConstantCache}, W_γdt)
nlcache.W_γdt = W_γdt
W_γdt
end
du_cache(nlsolver::AbstractNLSolver) = du_cache(nlsolver.cache)
du_cache(::AbstractNLSolverCache) = nothing
du_cache(nlcache::Union{NLFunctionalCache,NLAndersonCache,NLNewtonCache}) = (nlcache.k,)
function du_alias_or_new(nlsolver::AbstractNLSolver, rate_prototype)
_du_cache = du_cache(nlsolver)
if _du_cache === nothing
zero(rate_prototype)
else
first(_du_cache)
end
end
mutable struct DAEResidualJacobianWrapper{AD,F,pType,duType,uType,alphaType,gammaType,tmpType,uprevType,tType} <: Function
f::F
p::pType
tmp_du::duType
tmp_u::uType
α::alphaType
invγdt::gammaType
tmp::tmpType
uprev::uprevType
t::tType
function DAEResidualJacobianWrapper(alg,f,p,α,invγdt,tmp,uprev,t)
isautodiff = alg_autodiff(alg)
if isautodiff
tmp_du = DiffEqBase.dualcache(uprev)
tmp_u = DiffEqBase.dualcache(uprev)
else
tmp_du = similar(uprev)
tmp_u = similar(uprev)
end
new{isautodiff,typeof(f),typeof(p),typeof(tmp_du),typeof(tmp_u),typeof(α),typeof(invγdt),typeof(tmp),typeof(uprev),typeof(t)}(f,p,tmp_du,tmp_u,α,invγdt,tmp,uprev,t)
end
end
is_autodiff(m::DAEResidualJacobianWrapper{AD}) where AD = AD
function (m::DAEResidualJacobianWrapper)(out,x)
if is_autodiff(m)
tmp_du = DiffEqBase.get_tmp(m.tmp_du, x)
tmp_u = DiffEqBase.get_tmp(m.tmp_u, x)
else
tmp_du = m.tmp_du
tmp_u = m.tmp_u
end
@. tmp_du = (m.α * x + m.tmp) * m.invγdt
@. tmp_u = x + m.uprev
m.f(out, tmp_du, tmp_u, m.p, m.t)
end
mutable struct DAEResidualDerivativeWrapper{F,pType,alphaType,gammaType,tmpType,uprevType,tType} <: Function
f::F
p::pType
α::alphaType
invγdt::gammaType
tmp::tmpType
uprev::uprevType
t::tType
end
function (m::DAEResidualDerivativeWrapper)(x)
tmp_du = (m.α * x + m.tmp) * m.invγdt
tmp_u = x + m.uprev
m.f(tmp_du, tmp_u, m.p, m.t)
end
DiffEqBase.has_jac(f::DAEResidualJacobianWrapper) = DiffEqBase.has_jac(f.f)
DiffEqBase.has_Wfact(f::DAEResidualJacobianWrapper) = DiffEqBase.has_Wfact(f.f)
DiffEqBase.has_Wfact_t(f::DAEResidualJacobianWrapper) = DiffEqBase.has_Wfact_t(f.f)
DiffEqBase.has_jac(f::DAEResidualDerivativeWrapper) = DiffEqBase.has_jac(f.f)
DiffEqBase.has_Wfact(f::DAEResidualDerivativeWrapper) = DiffEqBase.has_Wfact(f.f)
DiffEqBase.has_Wfact_t(f::DAEResidualDerivativeWrapper) = DiffEqBase.has_Wfact_t(f.f)
function build_nlsolver(alg,u,uprev,p,t,dt,f,rate_prototype,uEltypeNoUnits,uBottomEltypeNoUnits,
tTypeNoUnits,γ,c,iip)
build_nlsolver(alg,u,uprev,p,t,dt,f,rate_prototype,uEltypeNoUnits,uBottomEltypeNoUnits,
tTypeNoUnits,γ,c,1,iip)
end
function build_nlsolver(alg,u,uprev,p,t,dt,f,rate_prototype,uEltypeNoUnits,uBottomEltypeNoUnits,
tTypeNoUnits,γ,c,α,iip)
build_nlsolver(alg,alg.nlsolve,u,uprev,p,t,dt,f,rate_prototype,uEltypeNoUnits,
uBottomEltypeNoUnits,tTypeNoUnits,γ,c,α,iip)
end
function build_nlsolver(alg,nlalg::Union{NLFunctional,NLAnderson,NLNewton},u,uprev,p,t,dt,
f,rate_prototype,uEltypeNoUnits,uBottomEltypeNoUnits,tTypeNoUnits,
γ,c,α,::Val{true})
#TODO
#nlalg = DiffEqBase.handle_defaults(alg, nlalg)
# define unitless type
uTolType = real(uBottomEltypeNoUnits)
isdae = alg isa DAEAlgorithm
# define fields of non-linear solver
z = zero(u); tmp = zero(u); ztmp = zero(u)
# build cache of non-linear solver
ustep = zero(u)
tstep = zero(t)
k = zero(rate_prototype)
atmp = similar(u, uEltypeNoUnits)
dz = zero(u)
if nlalg isa NLNewton
nf = nlsolve_f(f, alg)
if islinear(f)
du1 = rate_prototype
uf = nothing
jac_config = nothing
linsolve = alg.linsolve(Val{:init},nf,u)
else
du1 = zero(rate_prototype)
if isdae
uf = DAEResidualJacobianWrapper(alg,f,p,α,inv(γ*dt),tmp,uprev,t)
else
uf = build_uf(alg,nf,t,p,Val(true))
end
jac_config = build_jac_config(alg,nf,uf,du1,uprev,u,ztmp,dz)
linsolve = alg.linsolve(Val{:init},uf,u)
end
# TODO: check if the solver is iterative
weight = zero(u)
tType = typeof(t)
invγdt = inv(oneunit(t) * one(uTolType))
J, W = build_J_W(alg,u,uprev,p,t,dt,f,uEltypeNoUnits,Val(true))
nlcache = NLNewtonCache(ustep,tstep,k,atmp,dz,J,W,true,true,true,tType(dt),du1,uf,jac_config,
linsolve,weight,invγdt,tType(nlalg.new_W_dt_cutoff),t)
elseif nlalg isa NLFunctional
nlcache = NLFunctionalCache(ustep,tstep,k,atmp,dz)
elseif nlalg isa NLAnderson
max_history = min(nlalg.max_history, nlalg.max_iter, length(z))
Δz₊s = [zero(z) for i in 1:max_history]
Q = Matrix{uEltypeNoUnits}(undef, length(z), max_history)
R = Matrix{uEltypeNoUnits}(undef, max_history, max_history)
γs = Vector{uEltypeNoUnits}(undef, max_history)
dzold = zero(z)
z₊old = zero(z)
nlcache = NLAndersonCache(ustep,tstep,atmp,k,dz,dzold,z₊old,Δz₊s,Q,R,γs,0,
nlalg.aa_start,nlalg.droptol)
end
# build non-linear solver
ηold = one(t)
NLSolver{true,tTypeNoUnits}(
z,tmp,ztmp,γ,c,α,nlalg,nlalg.κ,
nlalg.fast_convergence_cutoff,ηold,0,nlalg.max_iter,Divergence,
nlcache)
end
function build_nlsolver(alg,nlalg::Union{NLFunctional,NLAnderson,NLNewton},u,uprev,p,t,dt,
f,rate_prototype,uEltypeNoUnits,uBottomEltypeNoUnits,tTypeNoUnits,
γ,c,α,::Val{false})
#TODO
#nlalg = DiffEqBase.handle_defaults(alg, nlalg)
# define unitless type
uTolType = real(uBottomEltypeNoUnits)
isdae = alg isa DAEAlgorithm
# define fields of non-linear solver
z = u; tmp = u; ztmp = u
# build cache of non-linear solver
tstep = zero(t)
if nlalg isa NLNewton
nf = nlsolve_f(f, alg)
if isdae
uf = DAEResidualDerivativeWrapper(f,p,α,inv(γ*dt),tmp,uprev,t)
else
uf = build_uf(alg,nf,t,p,Val(false))
end
tType = typeof(t)
invγdt = inv(oneunit(t) * one(uTolType))
J, W = build_J_W(alg,u,uprev,p,t,dt,f,uEltypeNoUnits,Val(false))
nlcache = NLNewtonConstantCache(tstep,J,W,true,true,true,tType(dt),uf,invγdt,tType(nlalg.new_W_dt_cutoff),t)
elseif nlalg isa NLFunctional
nlcache = NLFunctionalConstantCache(tstep)
elseif nlalg isa NLAnderson
max_history = min(nlalg.max_history, nlalg.max_iter, length(z))
Δz₊s = Vector{typeof(z)}(undef, max_history)
Q = Matrix{uEltypeNoUnits}(undef, length(z), max_history)
R = Matrix{uEltypeNoUnits}(undef, max_history, max_history)
γs = Vector{uEltypeNoUnits}(undef, max_history)
dz = u
dzold = u
z₊old = u
nlcache = NLAndersonConstantCache(tstep,dz,dzold,z₊old,Δz₊s,Q,R,γs,0,nlalg.aa_start,nlalg.droptol)
end
# build non-linear solver
ηold = one(tTypeNoUnits)
NLSolver{false,tTypeNoUnits}(
z, tmp, ztmp, γ, c, α, nlalg, nlalg.κ,
nlalg.fast_convergence_cutoff, ηold, 0, nlalg.max_iter, Divergence,
nlcache)
end
## Anderson acceleration
"""
anderson(z, cache)
Return the next iterate of the fixed-point iteration `z = g(z)` by performing Anderson
acceleration based on the current iterate `z` and the settings and history in the `cache`.
"""
@muladd function anderson(z, cache)
@unpack dz,Δz₊s,z₊old,dzold,R,Q,γs,history,droptol = cache
# increase size of history
history += 1
# remove oldest history if maximum size is exceeded
max_history = length(Δz₊s)
if history > max_history
# circularly shift differences of z₊
for i in 1:(max_history-1)
Δz₊s[i] = Δz₊s[i + 1]
end
# delete left-most column of QR decomposition
qrdelete!(Q, R, max_history)
# update size of history
history = max_history
end
# update history of differences of z₊
Δz₊s[history] = @.. z - z₊old
# replace/add difference of residuals as right-most column to QR decomposition
qradd!(Q, R, _vec(dz .- dzold), history)
# update cached values
cache.dzold = dz
cache.z₊old = z
# define current Q and R matrices
Qcur, Rcur = view(Q, :, 1:history), UpperTriangular(view(R, 1:history, 1:history))
# check condition (TODO: incremental estimation)
if droptol !== nothing
while cond(R) > droptol && history > 1
qrdelete!(Q, R, history)
history -= 1
Qcur, Rcur = view(Q, :, 1:history), UpperTriangular(view(R, 1:history, 1:history))
end
end
# save updated history
cache.history = history
# solve least squares problem
γscur = view(γs, 1:history)
ldiv!(Rcur, mul!(γscur, Qcur', _vec(dz)))
# update next iterate
for i in 1:history
z = @.. z - γs[i] * Δz₊s[i]
end
z
end
"""
anderson!(z, cache)
Update the current iterate `z` of the fixed-point iteration `z = g(z)` in-place
by performing Anderson acceleration based on the settings and history in the `cache`.
"""
@muladd function anderson!(z, cache)
@unpack dz,z₊old,dzold,Δz₊s,γs,R,Q,history,droptol = cache
# increase size of history
history += 1
# remove oldest history if maximum size is exceeded
max_history = length(Δz₊s)
if history > max_history
# circularly shift differences of z₊
ptr = Δz₊s[1]
for i in 1:(max_history-1)
Δz₊s[i] = Δz₊s[i + 1]
end
Δz₊s[max_history] = ptr
# delete left-most column of QR decomposition
qrdelete!(Q, R, max_history)
# update size of history
history = max_history
end
# update history of differences of z₊
@.. Δz₊s[history] = z - z₊old
# replace/add difference of residuals as right-most column to QR decomposition
@.. dzold = dz - dzold
qradd!(Q, R, vec(dzold), history)
# update cached values
@.. dzold = dz
@.. z₊old = z
# define current Q and R matrices
Qcur, Rcur = view(Q, :, 1:history), UpperTriangular(view(R, 1:history, 1:history))
# check condition (TODO: incremental estimation)
if droptol !== nothing
while cond(R) > droptol && history > 1
qrdelete!(Q, R, history)
history -= 1
Qcur, Rcur = view(Q, :, 1:history), UpperTriangular(view(R, 1:history, 1:history))
end
end
# save updated history
cache.history = history
# solve least squares problem
γscur = view(γs, 1:history)
ldiv!(Rcur, mul!(γscur, Qcur', vec(dz)))
# update next iterate
for i in 1:history
@.. z = z - γs[i] * Δz₊s[i]
end
nothing
end
## resize
function resize_nlsolver!(integrator::DiffEqBase.DEIntegrator, i::Int)
isdefined(integrator.cache, :nlsolver) || return
@unpack nlsolver = integrator.cache
if nlsolver isa AbstractArray
for idx in eachindex(nlsolver)
resize!(nlsolver[idx], integrator, i)
end
else
resize!(nlsolver, integrator, i)
end
nothing
end
function Base.resize!(nlsolver::AbstractNLSolver, integrator, i::Int)
resize!(nlsolver.z, i)
resize!(nlsolver.tmp, i)
resize!(nlsolver.ztmp, i)
resize!(nlsolver.cache, nlsolver, integrator, i)
end
## default: dispatch only on the cache
Base.resize!(cache::AbstractNLSolverCache, nlsolver, integrator, i::Int) =
Base.resize!(cache, i)