Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 2 additions & 0 deletions src/ModelPredictiveControl.jl
Original file line number Diff line number Diff line change
Expand Up @@ -9,6 +9,8 @@ using RecipesBase

using DifferentiationInterface: ADTypes.AbstractADType, AutoForwardDiff
using DifferentiationInterface: AutoSparse, SecondOrder
using DifferentiationInterface: gradient, jacobian, hessian
using DifferentiationInterface: value_and_gradient, value_gradient_and_hessian
using DifferentiationInterface: gradient!, value_and_gradient!, prepare_gradient
using DifferentiationInterface: jacobian!, value_and_jacobian!, prepare_jacobian
using DifferentiationInterface: hessian!, value_gradient_and_hessian!, prepare_hessian
Expand Down
13 changes: 11 additions & 2 deletions src/controller/execute.jl
Original file line number Diff line number Diff line change
Expand Up @@ -105,8 +105,17 @@ The function should be called after calling [`moveinput!`](@ref). It returns the
- `:d` : current measured disturbance, ``\mathbf{d}(k)``

For [`LinMPC`](@ref) and [`NonLinMPC`](@ref), the field `:sol` also contains the optimizer
solution summary that can be printed. Lastly, the economical cost `:JE` and the custom
nonlinear constraints `:gc` values at the optimum are also available for [`NonLinMPC`](@ref).
solution summary that can be printed. Lastly, for [`NonLinMPC`](@ref), the following fields
are also available:

- `:JE`: economic cost value at the optimum, ``J_E``
- `:gc`: custom nonlinear constraints values at the optimum, ``\mathbf{g_c}``
- `:∇J` or *`:nablaJ`* : gradient of the objective function, ``\mathbf{\nabla} J``
- `:∇²J` or *`:nabla2J`* : Hessian of the objective function, ``\mathbf{\nabla^2}J``
- `:∇g` or *`:nablag`* : Jacobian of the inequality constraint, ``\mathbf{\nabla g}``
- `:∇²ℓg` or *`:nabla2lg`* : Hessian of the inequality Lagrangian, ``\mathbf{\nabla^2}\ell_{\mathbf{g}}``
- `:∇geq` or *`:nablageq`* : Jacobian of the equality constraint, ``\mathbf{\nabla g_{eq}}``
- `:∇²ℓgeq` or *`:nabla2lgeq`* : Hessian of the equality Lagrangian, ``\mathbf{\nabla^2}\ell_{\mathbf{g_{eq}}}``

# Examples
```jldoctest
Expand Down
100 changes: 99 additions & 1 deletion src/controller/nonlinmpc.jl
Original file line number Diff line number Diff line change
Expand Up @@ -526,9 +526,11 @@ end
"""
addinfo!(info, mpc::NonLinMPC) -> info

For [`NonLinMPC`](@ref), add `:sol` and the optimal economic cost `:JE`.
For [`NonLinMPC`](@ref), add `:sol`, the custom nonlinear objective `:JE`, the custom
constraint `:gc`, and the various derivatives.
"""
function addinfo!(info, mpc::NonLinMPC{NT}) where NT<:Real
# --- variables specific to NonLinMPC ---
U, Ŷ, D̂, ŷ, d, ϵ = info[:U], info[:Ŷ], info[:D̂], info[:ŷ], info[:d], info[:ϵ]
Ue = [U; U[(end - mpc.estim.model.nu + 1):end]]
Ŷe = [ŷ; Ŷ]
Expand All @@ -539,6 +541,102 @@ function addinfo!(info, mpc::NonLinMPC{NT}) where NT<:Real
info[:JE] = JE
info[:gc] = LHS
info[:sol] = JuMP.solution_summary(mpc.optim, verbose=true)
# --- derivatives ---
model, optim = mpc.estim.model, mpc.optim
transcription = mpc.transcription
nu, ny, nx̂, nϵ = model.nu, model.ny, mpc.estim.nx̂, mpc.nϵ
nk = get_nk(model, transcription)
Hp, Hc = mpc.Hp, mpc.Hc
ng = length(mpc.con.i_g)
nc, neq = mpc.con.nc, mpc.con.neq
nU, nŶ, nX̂, nK = mpc.Hp*nu, Hp*ny, Hp*nx̂, Hp*nk
nΔŨ, nUe, nŶe = nu*Hc + nϵ, nU + nu, nŶ + ny
ΔŨ = zeros(NT, nΔŨ)
x̂0end = zeros(NT, nx̂)
K0 = zeros(NT, nK)
Ue, Ŷe = zeros(NT, nUe), zeros(NT, nŶe)
U0, Ŷ0 = zeros(NT, nU), zeros(NT, nŶ)
Û0, X̂0 = zeros(NT, nU), zeros(NT, nX̂)
gc, g = zeros(NT, nc), zeros(NT, ng)
geq = zeros(NT, neq)
J_cache = (
Cache(ΔŨ), Cache(x̂0end), Cache(Ue), Cache(Ŷe), Cache(U0), Cache(Ŷ0),
Cache(Û0), Cache(K0), Cache(X̂0),
Cache(gc), Cache(g), Cache(geq),
)
function J!(Z̃, ΔŨ, x̂0end, Ue, Ŷe, U0, Ŷ0, Û0, K0, X̂0, gc, g, geq)
update_predictions!(ΔŨ, x̂0end, Ue, Ŷe, U0, Ŷ0, Û0, K0, X̂0, gc, g, geq, mpc, Z̃)
return obj_nonlinprog!(Ŷ0, U0, mpc, model, Ue, Ŷe, ΔŨ)
end
if !isnothing(mpc.hessian)
_, ∇J, ∇²J = value_gradient_and_hessian(J!, mpc.hessian, mpc.Z̃, J_cache...)
else
∇J, ∇²J = gradient(J!, mpc.gradient, mpc.Z̃, J_cache...), nothing
end
JNT = typeof(optim).parameters[1]
nonlin_constraints = JuMP.all_constraints(
optim, JuMP.Vector{JuMP.VariableRef}, MOI.VectorNonlinearOracle{JNT}
)
g_con, geq_con = nonlin_constraints
display(g_con)
display(geq_con)
λ, λeq = JuMP.dual.(g_con), JuMP.dual.(geq_con)
println(JuMP.dual.(JuMP.FixBoundRef.(optim[:Z̃var])))
display(λ)
display(λeq)

∇g_cache = (
Cache(ΔŨ), Cache(x̂0end), Cache(Ue), Cache(Ŷe), Cache(U0), Cache(Ŷ0),
Cache(Û0), Cache(K0), Cache(X̂0),
Cache(gc), Cache(geq),
)
function g!(g, Z̃, ΔŨ, x̂0end, Ue, Ŷe, U0, Ŷ0, Û0, K0, X̂0, gc, geq)
update_predictions!(ΔŨ, x̂0end, Ue, Ŷe, U0, Ŷ0, Û0, K0, X̂0, gc, g, geq, mpc, Z̃)
return nothing
end
∇g = jacobian(g!, g, mpc.jacobian, mpc.Z̃, ∇g_cache...)
#=
if !isnothing(mpc.hessian)
function ℓ_g(Z̃, λ, ΔŨ, x̂0end, Ue, Ŷe, U0, Ŷ0, Û0, K0, X̂0, gc, geq, g)
update_predictions!(ΔŨ, x̂0end, Ue, Ŷe, U0, Ŷ0, Û0, K0, X̂0, gc, g, geq, mpc, Z̃)
return dot(λ, g)
end
∇²g_cache = (
Cache(ΔŨ), Cache(x̂0end), Cache(Ue), Cache(Ŷe), Cache(U0), Cache(Ŷ0),
Cache(Û0), Cache(K0), Cache(X̂0),
Cache(gc), Cache(geq), Cache(g)
)
∇²ℓg = hessian(ℓ_g, mpc.hessian, mpc.Z̃, Constant(λ), ∇²g_cache...)
else
∇²ℓg = nothing
end
=# ∇²ℓg = nothing #TODO: delete this line when enabling the above block
geq_cache = (
Cache(ΔŨ), Cache(x̂0end), Cache(Ue), Cache(Ŷe), Cache(U0), Cache(Ŷ0),
Cache(Û0), Cache(K0), Cache(X̂0),
Cache(gc), Cache(g)
)
function geq!(geq, Z̃, ΔŨ, x̂0end, Ue, Ŷe, U0, Ŷ0, Û0, K0, X̂0, gc, g)
update_predictions!(ΔŨ, x̂0end, Ue, Ŷe, U0, Ŷ0, Û0, K0, X̂0, gc, g, geq, mpc, Z̃)
return nothing
end
∇geq = jacobian(geq!, geq, mpc.jacobian, mpc.Z̃, geq_cache...)
∇²ℓgeq = nothing # TODO: implement later


info[:∇J] = ∇J
info[:∇²J] = ∇²J
info[:∇g] = ∇g
info[:∇²ℓg] = ∇²ℓg
info[:∇geq] = ∇geq
info[:∇²ℓgeq] = ∇²ℓgeq
# --- non-Unicode fields ---
info[:nablaJ] = ∇J
info[:nabla2J] = ∇²J
info[:nablag] = ∇g
info[:nabla2lg] = ∇²ℓg
info[:nablageq] = ∇geq
info[:nabla2lgeq] = ∇²ℓgeq
return info
end

Expand Down
4 changes: 1 addition & 3 deletions src/estimator/mhe/construct.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1504,9 +1504,7 @@ function get_nonlincon_oracle(
return dot(λi, gi)
end
Z̃_∇gi = fill(myNaN, nZ̃) # NaN to force update_predictions! at first call
∇gi_cache = (
Cache(V̂), Cache(X̂0), Cache(û0), Cache(k0), Cache(ŷ0), Cache(g)
)
∇gi_cache = (Cache(V̂), Cache(X̂0), Cache(û0), Cache(k0), Cache(ŷ0), Cache(g))
# temporarily "fill" the estimation window for the preparation of the gradient:
estim.Nk[] = He
∇gi_prep = prepare_jacobian(gi!, gi, jac, Z̃_∇gi, ∇gi_cache...; strict)
Expand Down
76 changes: 76 additions & 0 deletions src/estimator/mhe/execute.jl
Original file line number Diff line number Diff line change
Expand Up @@ -163,9 +163,85 @@ function getinfo(estim::MovingHorizonEstimator{NT}) where NT<:Real
info[:Yhatm] = info[:Ŷm]
# --- deprecated fields ---
info[:ϵ] = info[:ε]
info = addinfo!(info, estim, model)
return info
end


"""
addinfo!(info, estim::MovingHorizonEstimator, model::NonLinModel)

For [`NonLinModel`](@ref), add the various derivatives.
"""
function addinfo!(
info, estim::MovingHorizonEstimator{NT}, model::NonLinModel
) where NT <:Real
# --- derivatives ---
optim, con = estim.optim, estim.con
nx̂, nym, nŷ, nu, nk = estim.nx̂, estim.nym, model.ny, model.nu, model.nk
He = estim.He
ng = length(con.i_g)
nV̂, nX̂, ng = He*nym, He*nx̂, length(con.i_g)
V̂, X̂0 = zeros(NT, nV̂), zeros(NT, nX̂)
k0 = zeros(NT, nk)
û0, ŷ0 = zeros(NT, nu), zeros(NT, nŷ)
g = zeros(NT, ng)
x̄ = zeros(NT, nx̂)
J_cache = (
Cache(V̂), Cache(X̂0), Cache(û0), Cache(k0), Cache(ŷ0),
Cache(g),
Cache(x̄),
)
function J!(Z̃, V̂, X̂0, û0, k0, ŷ0, g, x̄)
update_prediction!(V̂, X̂0, û0, k0, ŷ0, g, estim, Z̃)
return obj_nonlinprog!(x̄, estim, model, V̂, Z̃)
end
if !isnothing(estim.hessian)
_, ∇J, ∇²J = value_gradient_and_hessian(J!, estim.hessian, estim.Z̃, J_cache...)
else
∇J, ∇²J = gradient(J!, estim.gradient, estim.Z̃, J_cache...), nothing
end
JNT = typeof(optim).parameters[1]
nonlin_constraints = JuMP.all_constraints(
optim, JuMP.Vector{JuMP.VariableRef}, MOI.VectorNonlinearOracle{JNT}
)
g_con = nonlin_constraints[1]
λ = JuMP.dual.(g_con)
display(λ)
∇g_cache = (Cache(V̂), Cache(X̂0), Cache(û0), Cache(k0), Cache(ŷ0))
function g!(g, Z̃, V̂, X̂0, û0, k0, ŷ0)
update_prediction!(V̂, X̂0, û0, k0, ŷ0, g, estim, Z̃)
return nothing
end
∇g = jacobian(g!, g, estim.jacobian, estim.Z̃, ∇g_cache...)
#=
if !isnothing(estim.hessian)
function ℓ_g(Z̃, λ, V̂, X̂0, û0, k0, ŷ0, g)
update_prediction!(V̂, X̂0, û0, k0, ŷ0, g, estim, Z̃)
return dot(λ, g)
end
∇²g_cache = (Cache(V̂), Cache(X̂0), Cache(û0), Cache(k0), Cache(ŷ0), Cache(g))
∇²ℓg = hessian(ℓ_g, estim.hessian, estim.Z̃, Constant(λ), ∇²g_cache...)
else
∇²ℓg = nothing
end
=# ∇²ℓg = nothing #TODO: delete this line when enabling the above block

info[:∇J] = ∇J
info[:∇²J] = ∇²J
info[:∇g] = ∇g
info[:∇²ℓg] = ∇²ℓg
# --- non-Unicode fields ---
info[:nablaJ] = ∇J
info[:nabla2J] = ∇²J
info[:nablag] = ∇g
info[:nabla2lg] = ∇²ℓg
return info
end

"Nothing to add in the `info` dict for [`LinModel`](@ref)."
addinfo!(info, ::MovingHorizonEstimator, ::LinModel) = info

"""
getε(estim::MovingHorizonEstimator, Z̃) -> ε

Expand Down
17 changes: 15 additions & 2 deletions src/general.jl
Original file line number Diff line number Diff line change
Expand Up @@ -16,6 +16,14 @@ const ALL_COLORING_ORDERS = (
RandomOrder(StableRNG(0), 0)
)

const HIDDEN_GETINFO_KEYS_MHE = (
:What, :xhatarr, :epsilon, :Xhat, :xhat, :Vhat, :Pbar, :xbar, :Yhat, :Yhatm, :ϵ
)

const HIDDEN_GETINFO_KEYS_MPC = (
:DeltaU, :epsilon, :Dhat, :yhat, :Yhat, :xhatend, :Yhats, :Rhaty, :Rhatu
)

"Termination status that means 'no solution available'."
const ERROR_STATUSES = (
JuMP.INFEASIBLE, JuMP.DUAL_INFEASIBLE, JuMP.LOCALLY_INFEASIBLE,
Expand All @@ -40,12 +48,17 @@ function info2debugstr(info)
mystr = "Content of getinfo dictionary:\n"
for (key, value) in info
(key == :sol) && continue
if key in HIDDEN_GETINFO_KEYS_MHE || key in HIDDEN_GETINFO_KEYS_MPC
# skip the redundant non-Unicode keys
continue
end
mystr *= " :$key => $value\n"
end
if haskey(info, :sol)
split_sol = split(string(info[:sol]), "\n")
solstr = join((lpad(line, length(line) + 2) for line in split_sol), "\n", "")
mystr *= " :sol => \n"*solstr
# Add the treeview prefix to each line
solstr = join((" " * line for line in split_sol), "\n")
mystr *= " :sol => \n" * solstr * "\n" # Ensure a trailing newline
end
return mystr
end
Expand Down
Loading