|
| 1 | +""" |
| 2 | + LimitedMemoryBroyden(; max_resets::Int = 3, linesearch = LineSearch(), |
| 3 | + threshold::Int = 10, reset_tolerance = nothing) |
| 4 | +
|
| 5 | +An implementation of `LimitedMemoryBroyden` with reseting and line search. |
| 6 | +
|
| 7 | +## Arguments |
| 8 | +
|
| 9 | + - `max_resets`: the maximum number of resets to perform. Defaults to `3`. |
| 10 | + - `reset_tolerance`: the tolerance for the reset check. Defaults to |
| 11 | + `sqrt(eps(eltype(u)))`. |
| 12 | + - `threshold`: the number of vectors to store in the low rank approximation. Defaults |
| 13 | + to `10`. |
| 14 | + - `linesearch`: the line search algorithm to use. Defaults to [`LineSearch()`](@ref), |
| 15 | + which means that no line search is performed. Algorithms from `LineSearches.jl` can be |
| 16 | + used here directly, and they will be converted to the correct `LineSearch`. It is |
| 17 | + recommended to use [LiFukushimaLineSearchCache](@ref) -- a derivative free linesearch |
| 18 | + specifically designed for Broyden's method. |
| 19 | +""" |
| 20 | +@concrete struct LimitedMemoryBroyden <: AbstractNewtonAlgorithm{false, Nothing} |
| 21 | + max_resets::Int |
| 22 | + threshold::Int |
| 23 | + linesearch |
| 24 | + reset_tolerance |
| 25 | +end |
| 26 | + |
| 27 | +function LimitedMemoryBroyden(; max_resets::Int = 3, linesearch = LineSearch(), |
| 28 | + threshold::Int = 10, reset_tolerance = nothing) |
| 29 | + linesearch = linesearch isa LineSearch ? linesearch : LineSearch(; method = linesearch) |
| 30 | + return LimitedMemoryBroyden(max_resets, threshold, linesearch, reset_tolerance) |
| 31 | +end |
| 32 | + |
| 33 | +@concrete mutable struct LimitedMemoryBroydenCache{iip} <: AbstractNonlinearSolveCache{iip} |
| 34 | + f |
| 35 | + alg |
| 36 | + u |
| 37 | + du |
| 38 | + fu |
| 39 | + fu2 |
| 40 | + dfu |
| 41 | + p |
| 42 | + U |
| 43 | + Vᵀ |
| 44 | + Ux |
| 45 | + xᵀVᵀ |
| 46 | + u_cache |
| 47 | + vᵀ_cache |
| 48 | + force_stop::Bool |
| 49 | + resets::Int |
| 50 | + iterations_since_reset::Int |
| 51 | + max_resets::Int |
| 52 | + maxiters::Int |
| 53 | + internalnorm |
| 54 | + retcode::ReturnCode.T |
| 55 | + abstol |
| 56 | + reset_tolerance |
| 57 | + reset_check |
| 58 | + prob |
| 59 | + stats::NLStats |
| 60 | + lscache |
| 61 | +end |
| 62 | + |
| 63 | +get_fu(cache::LimitedMemoryBroydenCache) = cache.fu |
| 64 | + |
| 65 | +function SciMLBase.__init(prob::NonlinearProblem{uType, iip}, alg::LimitedMemoryBroyden, |
| 66 | + args...; alias_u0 = false, maxiters = 1000, abstol = 1e-6, internalnorm = DEFAULT_NORM, |
| 67 | + kwargs...) where {uType, iip} |
| 68 | + @unpack f, u0, p = prob |
| 69 | + u = alias_u0 ? u0 : deepcopy(u0) |
| 70 | + if u isa Number |
| 71 | + # If u is a number then we simply use Broyden |
| 72 | + return SciMLBase.__init(prob, |
| 73 | + GeneralBroyden(; alg.max_resets, alg.reset_tolerance, |
| 74 | + alg.linesearch), args...; alias_u0, maxiters, abstol, internalnorm, kwargs...) |
| 75 | + end |
| 76 | + fu = evaluate_f(prob, u) |
| 77 | + threshold = min(alg.threshold, maxiters) |
| 78 | + U, Vᵀ = __init_low_rank_jacobian(u, fu, threshold) |
| 79 | + du = -fu |
| 80 | + reset_tolerance = alg.reset_tolerance === nothing ? sqrt(eps(eltype(u))) : |
| 81 | + alg.reset_tolerance |
| 82 | + reset_check = x -> abs(x) ≤ reset_tolerance |
| 83 | + return LimitedMemoryBroydenCache{iip}(f, alg, u, du, fu, zero(fu), |
| 84 | + zero(fu), p, U, Vᵀ, similar(u, threshold), similar(u, 1, threshold), |
| 85 | + zero(u), zero(u), false, 0, 0, alg.max_resets, maxiters, internalnorm, |
| 86 | + ReturnCode.Default, abstol, reset_tolerance, reset_check, prob, |
| 87 | + NLStats(1, 0, 0, 0, 0), |
| 88 | + init_linesearch_cache(alg.linesearch, f, u, p, fu, Val(iip))) |
| 89 | +end |
| 90 | + |
| 91 | +function perform_step!(cache::LimitedMemoryBroydenCache{true}) |
| 92 | + @unpack f, p, du, u = cache |
| 93 | + T = eltype(u) |
| 94 | + |
| 95 | + α = perform_linesearch!(cache.lscache, u, du) |
| 96 | + axpy!(α, du, u) |
| 97 | + f(cache.fu2, u, p) |
| 98 | + |
| 99 | + cache.internalnorm(cache.fu2) < cache.abstol && (cache.force_stop = true) |
| 100 | + cache.stats.nf += 1 |
| 101 | + |
| 102 | + cache.force_stop && return nothing |
| 103 | + |
| 104 | + # Update the Inverse Jacobian Approximation |
| 105 | + cache.dfu .= cache.fu2 .- cache.fu |
| 106 | + |
| 107 | + # Only try to reset if we have enough iterations since last reset |
| 108 | + if cache.iterations_since_reset > size(cache.U, 1) && |
| 109 | + (all(cache.reset_check, du) || all(cache.reset_check, cache.dfu)) |
| 110 | + if cache.resets ≥ cache.max_resets |
| 111 | + cache.retcode = ReturnCode.Unstable |
| 112 | + cache.force_stop = true |
| 113 | + return nothing |
| 114 | + end |
| 115 | + cache.iterations_since_reset = 0 |
| 116 | + cache.resets += 1 |
| 117 | + cache.du .= -cache.fu |
| 118 | + else |
| 119 | + idx = min(cache.iterations_since_reset, size(cache.U, 1)) |
| 120 | + U_part = selectdim(cache.U, 1, 1:idx) |
| 121 | + Vᵀ_part = selectdim(cache.Vᵀ, 2, 1:idx) |
| 122 | + |
| 123 | + __lbroyden_matvec!(_vec(cache.vᵀ_cache), cache.Ux, U_part, Vᵀ_part, _vec(cache.du)) |
| 124 | + __lbroyden_rmatvec!(_vec(cache.u_cache), cache.xᵀVᵀ, U_part, Vᵀ_part, |
| 125 | + _vec(cache.dfu)) |
| 126 | + cache.u_cache .= (du .- cache.u_cache) ./ |
| 127 | + (dot(cache.vᵀ_cache, cache.dfu) .+ T(1e-5)) |
| 128 | + |
| 129 | + idx = mod1(cache.iterations_since_reset + 1, size(cache.U, 1)) |
| 130 | + selectdim(cache.U, 1, idx) .= _vec(cache.u_cache) |
| 131 | + selectdim(cache.Vᵀ, 2, idx) .= _vec(cache.vᵀ_cache) |
| 132 | + |
| 133 | + idx = min(cache.iterations_since_reset + 1, size(cache.U, 1)) |
| 134 | + U_part = selectdim(cache.U, 1, 1:idx) |
| 135 | + Vᵀ_part = selectdim(cache.Vᵀ, 2, 1:idx) |
| 136 | + __lbroyden_matvec!(_vec(cache.du), cache.Ux, U_part, Vᵀ_part, _vec(cache.fu2)) |
| 137 | + cache.du .*= -1 |
| 138 | + cache.iterations_since_reset += 1 |
| 139 | + end |
| 140 | + |
| 141 | + cache.fu .= cache.fu2 |
| 142 | + |
| 143 | + return nothing |
| 144 | +end |
| 145 | + |
| 146 | +function perform_step!(cache::LimitedMemoryBroydenCache{false}) |
| 147 | + @unpack f, p = cache |
| 148 | + T = eltype(cache.u) |
| 149 | + |
| 150 | + α = perform_linesearch!(cache.lscache, cache.u, cache.du) |
| 151 | + cache.u = cache.u .+ α * cache.du |
| 152 | + cache.fu2 = f(cache.u, p) |
| 153 | + |
| 154 | + cache.internalnorm(cache.fu2) < cache.abstol && (cache.force_stop = true) |
| 155 | + cache.stats.nf += 1 |
| 156 | + |
| 157 | + cache.force_stop && return nothing |
| 158 | + |
| 159 | + # Update the Inverse Jacobian Approximation |
| 160 | + cache.dfu .= cache.fu2 .- cache.fu |
| 161 | + |
| 162 | + # Only try to reset if we have enough iterations since last reset |
| 163 | + if cache.iterations_since_reset > size(cache.U, 1) && |
| 164 | + (all(cache.reset_check, cache.du) || all(cache.reset_check, cache.dfu)) |
| 165 | + if cache.resets ≥ cache.max_resets |
| 166 | + cache.retcode = ReturnCode.Unstable |
| 167 | + cache.force_stop = true |
| 168 | + return nothing |
| 169 | + end |
| 170 | + cache.iterations_since_reset = 0 |
| 171 | + cache.resets += 1 |
| 172 | + cache.du = -cache.fu |
| 173 | + else |
| 174 | + idx = min(cache.iterations_since_reset, size(cache.U, 1)) |
| 175 | + U_part = selectdim(cache.U, 1, 1:idx) |
| 176 | + Vᵀ_part = selectdim(cache.Vᵀ, 2, 1:idx) |
| 177 | + |
| 178 | + cache.vᵀ_cache = _restructure(cache.vᵀ_cache, |
| 179 | + __lbroyden_matvec(U_part, Vᵀ_part, _vec(cache.du))) |
| 180 | + cache.u_cache = _restructure(cache.u_cache, |
| 181 | + __lbroyden_rmatvec(U_part, Vᵀ_part, _vec(cache.dfu))) |
| 182 | + cache.u_cache = (cache.du .- cache.u_cache) ./ |
| 183 | + (dot(cache.vᵀ_cache, cache.dfu) .+ T(1e-5)) |
| 184 | + |
| 185 | + idx = mod1(cache.iterations_since_reset + 1, size(cache.U, 1)) |
| 186 | + selectdim(cache.U, 1, idx) .= _vec(cache.u_cache) |
| 187 | + selectdim(cache.Vᵀ, 2, idx) .= _vec(cache.vᵀ_cache) |
| 188 | + |
| 189 | + idx = min(cache.iterations_since_reset + 1, size(cache.U, 1)) |
| 190 | + U_part = selectdim(cache.U, 1, 1:idx) |
| 191 | + Vᵀ_part = selectdim(cache.Vᵀ, 2, 1:idx) |
| 192 | + cache.du = _restructure(cache.du, |
| 193 | + -__lbroyden_matvec(U_part, Vᵀ_part, _vec(cache.fu2))) |
| 194 | + cache.iterations_since_reset += 1 |
| 195 | + end |
| 196 | + |
| 197 | + cache.fu = cache.fu2 |
| 198 | + |
| 199 | + return nothing |
| 200 | +end |
| 201 | + |
| 202 | +function SciMLBase.reinit!(cache::LimitedMemoryBroydenCache{iip}, u0 = cache.u; p = cache.p, |
| 203 | + abstol = cache.abstol, maxiters = cache.maxiters) where {iip} |
| 204 | + cache.p = p |
| 205 | + if iip |
| 206 | + recursivecopy!(cache.u, u0) |
| 207 | + cache.f(cache.fu, cache.u, p) |
| 208 | + else |
| 209 | + # don't have alias_u0 but cache.u is never mutated for OOP problems so it doesn't matter |
| 210 | + cache.u = u0 |
| 211 | + cache.fu = cache.f(cache.u, p) |
| 212 | + end |
| 213 | + cache.abstol = abstol |
| 214 | + cache.maxiters = maxiters |
| 215 | + cache.stats.nf = 1 |
| 216 | + cache.stats.nsteps = 1 |
| 217 | + cache.resets = 0 |
| 218 | + cache.iterations_since_reset = 0 |
| 219 | + cache.force_stop = false |
| 220 | + cache.retcode = ReturnCode.Default |
| 221 | + return cache |
| 222 | +end |
| 223 | + |
| 224 | +@views function __lbroyden_matvec!(y::AbstractVector, Ux::AbstractVector, |
| 225 | + U::AbstractMatrix, Vᵀ::AbstractMatrix, x::AbstractVector) |
| 226 | + # Computes Vᵀ × U × x |
| 227 | + η = size(U, 1) |
| 228 | + if η == 0 |
| 229 | + y .= x |
| 230 | + return nothing |
| 231 | + end |
| 232 | + mul!(Ux[1:η], U, x) |
| 233 | + mul!(y, Vᵀ[:, 1:η], Ux[1:η]) |
| 234 | + return nothing |
| 235 | +end |
| 236 | + |
| 237 | +@views function __lbroyden_matvec(U::AbstractMatrix, Vᵀ::AbstractMatrix, x::AbstractVector) |
| 238 | + # Computes Vᵀ × U × x |
| 239 | + size(U, 1) == 0 && return x |
| 240 | + return Vᵀ * (U * x) |
| 241 | +end |
| 242 | + |
| 243 | +@views function __lbroyden_rmatvec!(y::AbstractVector, xᵀVᵀ::AbstractMatrix, |
| 244 | + U::AbstractMatrix, Vᵀ::AbstractMatrix, x::AbstractVector) |
| 245 | + # Computes xᵀ × Vᵀ × U |
| 246 | + η = size(U, 1) |
| 247 | + if η == 0 |
| 248 | + y .= x |
| 249 | + return nothing |
| 250 | + end |
| 251 | + mul!(xᵀVᵀ[:, 1:η], x', Vᵀ) |
| 252 | + mul!(y', xᵀVᵀ[:, 1:η], U) |
| 253 | + return nothing |
| 254 | +end |
| 255 | + |
| 256 | +@views function __lbroyden_rmatvec(U::AbstractMatrix, Vᵀ::AbstractMatrix, x::AbstractVector) |
| 257 | + # Computes xᵀ × Vᵀ × U |
| 258 | + size(U, 1) == 0 && return x |
| 259 | + return (x' * Vᵀ) * U |
| 260 | +end |
0 commit comments