Manual: Nonlinear Design

Nonlinear Model

In this example, the goal is to control the angular position $θ$ of a pendulum attached to a motor. Knowing that the manipulated input is the motor torque $τ$, the I/O vectors are:

\[\begin{aligned} \mathbf{u} &= τ \\ \mathbf{y} &= θ \end{aligned}\]

The following figure presents the system:

pendulum

The plant model is nonlinear:

\[\begin{aligned} \dot{θ}(t) &= ω(t) \\ \dot{ω}(t) &= -\frac{g}{L}\sin\big( θ(t) \big) - \frac{K}{m} ω(t) + \frac{1}{m L^2} τ(t) \end{aligned}\]

in which $g$ is the gravitational acceleration in m/s², $L$, the pendulum length in m, $K$, the friction coefficient at the pivot point in kg/s, and $m$, the mass attached at the end of the pendulum in kg. The NonLinModel constructor assumes by default that the state function f is continuous in time, that is, an ordinary differential equation system (like here):

using ModelPredictiveControl
function pendulum(par, x, u)
    g, L, K, m = par        # [m/s²], [m], [kg/s], [kg]
    θ, ω = x[1], x[2]       # [rad], [rad/s]
    τ  = u[1]               # [N m]
    dθ = ω
    dω = -g/L*sin(θ) - K/m*ω + τ/m/L^2
    return [dθ, dω]
end
# declared constants, to avoid type-instability in the f function, for speed:
const par = (9.8, 0.4, 1.2, 0.3)
f(x, u, _ ) = pendulum(par, x, u)
h(x, _ )    = [180/π*x[1]]  # [°]
Ts, nu, nx, ny = 0.1, 1, 2, 1
model = NonLinModel(f, h, Ts, nu, nx, ny)
NonLinModel with a sample time Ts = 0.1 s, RungeKutta solver and:
 1 manipulated inputs u
 2 states x
 1 outputs y
 0 measured disturbances d

The output function $\mathbf{h}$ converts the $θ$ angle to degrees. Note that special characters like $θ$ can be typed in the Julia REPL or VS Code by typing \theta and pressing the <TAB> key. The tuple par is constant here to improve the performance. Note that a 4th order RungeKutta differential equation solver is used by default. It is good practice to first simulate model using sim! as a quick sanity check:

using Plots
u = [0.5]
N = 35
res = sim!(model, N, u)
plot(res, plotu=false)

plot1_NonLinMPC

Nonlinear Model Predictive Controller

An UnscentedKalmanFilter estimates the plant state :

α=0.01; σQ=[0.1, 0.5]; σR=[0.5]; nint_u=[1]; σQint_u=[0.1]
estim = UnscentedKalmanFilter(model; α, σQ, σR, nint_u, σQint_u)
UnscentedKalmanFilter estimator with a sample time Ts = 0.1 s, NonLinModel and:
 1 manipulated inputs u (1 integrating states)
 3 estimated states x̂
 1 measured outputs ym (0 integrating states)
 0 unmeasured outputs yu
 0 measured disturbances d

The vectors σQ and σR σR are the standard deviations of the process and sensor noises, respectively. The value for the velocity $ω$ is higher here (σQ second value) since $\dot{ω}(t)$ equation includes an uncertain parameter: the friction coefficient $K$. Also, the argument nint_u explicitly adds one integrating state at the model input, the motor torque $τ$, with an associated standard deviation σQint_u of 0.1 N m. The estimator tuning is tested on a plant with a 25 % larger friction coefficient $K$:

const par_plant = (par[1], par[2], 1.25*par[3], par[4])
f_plant(x, u, _ ) = pendulum(par_plant, x, u)
plant = NonLinModel(f_plant, h, Ts, nu, nx, ny)
res = sim!(estim, N, [0.5], plant=plant, y_noise=[0.5])
plot(res, plotu=false, plotxwithx̂=true)

plot2_NonLinMPC

The estimate $x̂_3$ is the integrating state on the torque $τ$ that compensates for static errors. The Kalman filter performance seems sufficient for control.

As the motor torque is limited to -1.5 to 1.5 N m, we incorporate the input constraints in a NonLinMPC:

Hp, Hc, Mwt, Nwt = 20, 2, [0.5], [2.5]
nmpc = NonLinMPC(estim; Hp, Hc, Mwt, Nwt, Cwt=Inf)
umin, umax = [-1.5], [+1.5]
nmpc = setconstraint!(nmpc; umin, umax)
NonLinMPC controller with a sample time Ts = 0.1 s, Ipopt optimizer, UnscentedKalmanFilter estimator and:
 20 prediction steps Hp
  2 control steps Hc
  0 slack variable ϵ (control constraints)
  1 manipulated inputs u (1 integrating states)
  3 estimated states x̂
  1 measured outputs ym (0 integrating states)
  0 unmeasured outputs yu
  0 measured disturbances d

The option Cwt=Inf disables the slack variable ϵ for constraint softening. We test mpc performance on plant by imposing an angular setpoint of 180° (inverted position):

res_ry = sim!(nmpc, N, [180.0], plant=plant, x_0=[0, 0], x̂_0=[0, 0, 0])
plot(res_ry)

plot3_NonLinMPC

The controller seems robust enough to variations on $K$ coefficient. Starting from this inverted position, the closed-loop response to a step disturbances of 10° is also satisfactory:

res_yd = sim!(nmpc, N, [180.0], plant=plant, x_0=[π, 0], x̂_0=[π, 0, 0], y_step=[10])
plot(res_yd)

plot4_NonLinMPC

Economic Model Predictive Controller

Economic MPC can achieve the same objective but with lower economical costs. For this case study, the controller will aim to reduce the energy consumed by the motor. The power (W) transmitted by the motor to the pendulum is:

\[Ẇ(t) = τ(t) ω(t)\]

Thus, the work (J) done by the motor from $t = t_0$ to $t_{end}$ is:

\[W = \int_{t_0}^{t_{end}} Ẇ(t) \mathrm{d}t = \int_{t_0}^{t_{end}} τ(t) ω(t) \mathrm{d}t\]

With the sampling time $T_s$ in s, the prediction horizon $H_p$, the limits defined as $t_0 = k T_s$ and $t_{end} = (k+H_p) T_s$, and the left-endpoint rectangle method for the integral, we get:

\[W ≈ T_s \sum_{j=0}^{H_p-1} τ(k + j) ω(k + j)\]

The objective function will now include an additive term that penalizes the work done by the motor $W$ to reduce the energy consumption. Notice that $W$ is a function of the manipulated input $τ$ and the angular speed $ω$, a state that is not measured (only the angle $θ$ is measured here). As the arguments of NonLinMPC economic function JE do not include the states, the speed is now defined as an unmeasured output to design a Kalman Filter similar to the previous one ($\mathbf{y^m} = θ$ and $\mathbf{y^u} = ω$):

h2(x, _ ) = [180/π*x[1], x[2]]
nu, nx, ny = 1, 2, 2
model2 = NonLinModel(f      , h2, Ts, nu, nx, ny)
plant2 = NonLinModel(f_plant, h2, Ts, nu, nx, ny)
estim2 = UnscentedKalmanFilter(model2; σQ, σR, nint_u, σQint_u, i_ym=[1])
UnscentedKalmanFilter estimator with a sample time Ts = 0.1 s, NonLinModel and:
 1 manipulated inputs u (1 integrating states)
 3 estimated states x̂
 1 measured outputs ym (0 integrating states)
 1 unmeasured outputs yu
 0 measured disturbances d

The plant2 object based on h2 is also required since sim! expects that the output vector of plant argument corresponds to the model output vector in mpc argument. We can now define the $J_E$ function and the empc controller:

function JE(UE, ŶE, _ )
    τ, ω = UE[1:end-1], ŶE[2:2:end-1]
    return Ts*sum(τ.*ω)
end
empc = NonLinMPC(estim2; Hp, Hc, Nwt, Mwt=[0.5, 0], Cwt=Inf, Ewt=3.5e3, JE=JE)
empc = setconstraint!(empc; umin, umax)
NonLinMPC controller with a sample time Ts = 0.1 s, Ipopt optimizer, UnscentedKalmanFilter estimator and:
 20 prediction steps Hp
  2 control steps Hc
  0 slack variable ϵ (control constraints)
  1 manipulated inputs u (1 integrating states)
  3 estimated states x̂
  1 measured outputs ym (0 integrating states)
  1 unmeasured outputs yu
  0 measured disturbances d

The keyword argument Ewt weights the economic costs relative to the other terms in the objective function. The term must be large enough to be significant but a too high value can lead to a static error on the angle setpoint. The second element of Mwt is zero since the speed $ω$ is not requested to track a setpoint. The closed-loop response to a 180° setpoint is similar:

res2_ry = sim!(empc, N, [180, 0], plant=plant2, x_0=[0, 0], x̂_0=[0, 0, 0])
plot(res2_ry)

plot5_NonLinMPC

and the energy consumption is slightly lower:

function calcW(res)
    τ, ω = res.U_data[1, 1:end-1], res.X_data[2, 1:end-1]
    return Ts*sum(τ.*ω)
end
Dict(:W_nmpc => calcW(res_ry), :W_empc => calcW(res2_ry))
Dict{Symbol, Float64} with 2 entries:
  :W_empc => 3.89437
  :W_nmpc => 3.91597

Also, for a 10° step disturbance:

res2_yd = sim!(empc, N, [180; 0]; plant=plant2, x_0=[π, 0], x̂_0=[π, 0, 0], y_step=[10, 0])
plot(res2_yd)

plot6_NonLinMPC

the new controller is able to recuperate more energy from the pendulum (i.e. negative work):

Dict(:W_nmpc => calcW(res_yd), :W_empc => calcW(res2_yd))
Dict{Symbol, Float64} with 2 entries:
  :W_empc => -0.149044
  :W_nmpc => -0.112247

Of course, this gain is only exploitable if the motor electronic includes some kind of regenerative circuitry.

Model Linearization

Nonlinear MPC is more computationally expensive than LinMPC. Solving the problem should always be faster than the sampling time $T_s = 0.1$ s for real-time operation. This requirement is sometimes hard to meet on electronics or mechanical systems because of the fast dynamics. To ease the design and comparison with LinMPC, the linearize function allows automatic linearization of NonLinModel based on ForwardDiff.jl. We first linearize model at the point $θ = π$ rad and $ω = τ = 0$ (inverted position):

linmodel = linearize(model, x=[π, 0], u=[0])
LinModel with a sample time Ts = 0.1 s and:
 1 manipulated inputs u
 2 states x
 1 outputs y
 0 measured disturbances d

A SteadyKalmanFilter and a LinMPC are designed from linmodel:

skf = SteadyKalmanFilter(linmodel; σQ, σR, nint_u, σQint_u)
mpc = LinMPC(skf; Hp, Hc, Mwt, Nwt, Cwt=Inf)
mpc = setconstraint!(mpc, umin=[-1.5], umax=[+1.5])
LinMPC controller with a sample time Ts = 0.1 s, OSQP optimizer, SteadyKalmanFilter estimator and:
 20 prediction steps Hp
  2 control steps Hc
  0 slack variable ϵ (control constraints)
  1 manipulated inputs u (1 integrating states)
  3 estimated states x̂
  1 measured outputs ym (0 integrating states)
  0 unmeasured outputs yu
  0 measured disturbances d

The linear controller has difficulties to reject the 10° step disturbance:

res_lin = sim!(mpc, N, [180.0]; plant, x_0=[π, 0], y_step=[10])
plot(res_lin)

plot7_NonLinMPC

Solving the optimization problem of mpc with DAQP optimizer instead of the default OSQP solver can help here. It is indeed documented that DAQP can perform better on small/medium dense matrices and unstable poles[1], which is obviously the case here (absolute value of unstable poles are greater than one):

using LinearAlgebra; poles = eigvals(linmodel.A)
2-element Vector{Float64}:
 0.48163285653799215
 1.3963025601286745

To install the solver, run:

using Pkg; Pkg.add("DAQP")

Constructing a LinMPC with DAQP:

using JuMP, DAQP
daqp = Model(DAQP.Optimizer, add_bridges=false)
mpc2 = LinMPC(skf; Hp, Hc, Mwt, Nwt, Cwt=Inf, optim=daqp)
mpc2 = setconstraint!(mpc2; umin, umax)
LinMPC controller with a sample time Ts = 0.1 s, DAQP optimizer, SteadyKalmanFilter estimator and:
 20 prediction steps Hp
  2 control steps Hc
  0 slack variable ϵ (control constraints)
  1 manipulated inputs u (1 integrating states)
  3 estimated states x̂
  1 measured outputs ym (0 integrating states)
  0 unmeasured outputs yu
  0 measured disturbances d

does improve the rejection of the step disturbance:

res_lin2 = sim!(mpc2, N, [180.0]; plant, x_0=[π, 0], y_step=[10])
plot(res_lin2)

plot8_NonLinMPC

The closed-loop performance is still lower than the nonlinear controller, as expected, but computations are about 210 times faster (0.000071 s versus 0.015 s per time steps, on average). However, remember that linmodel is only valid for angular positions near 180°. For example, the 180° setpoint response from 0° is unsatisfactory since the predictions are poor in the first quadrant:

res_lin3 = sim!(mpc2, N, [180.0]; plant, x_0=[0, 0])
plot(res_lin3)

plot9_NonLinMPC

Multiple linearized model and controller objects are required for large deviations from this operating point. This is known as gain scheduling. Another approach is adapting the model of the LinMPC instance based on repeated online linearization.

Adapting the Model via Successive Linearization

The setmodel! method allows online adaptation of a linear plant model. Combined with the automatic linearization of linearize, a successive linearization MPC can be designed with minimal efforts. The SteadyKalmanFilter does not support setmodel!, so we need to use the time-varying KalmanFilter instead:

kf   = KalmanFilter(linmodel; σQ, σR, nint_u, σQint_u)
mpc3 = LinMPC(kf; Hp, Hc, Mwt, Nwt, Cwt=Inf, optim=daqp)
mpc3 = setconstraint!(mpc3; umin, umax)
LinMPC controller with a sample time Ts = 0.1 s, DAQP optimizer, KalmanFilter estimator and:
 20 prediction steps Hp
  2 control steps Hc
  0 slack variable ϵ (control constraints)
  1 manipulated inputs u (1 integrating states)
  3 estimated states x̂
  1 measured outputs ym (0 integrating states)
  0 unmeasured outputs yu
  0 measured disturbances d

We create a function that simulates the plant and the adaptive controller:

function test_slmpc(nonlinmodel, mpc, ry, plant; x_0=plant.xop, y_step=0)
    N = 35
    U_data, Y_data, Ry_data = zeros(plant.nu, N), zeros(plant.ny, N), zeros(plant.ny, N)
    setstate!(plant, x_0)
    u, y = [0.0], plant()
    x̂ = initstate!(mpc, u, y)
    for i = 1:N
        y = plant() .+ y_step
        u = moveinput!(mpc, ry)
        linmodel = linearize(nonlinmodel; u, x=x̂[1:2])
        setmodel!(mpc, linmodel)
        U_data[:,i], Y_data[:,i], Ry_data[:,i] = u, y, ry
        x̂ = updatestate!(mpc, u, y) # update mpc state estimate
        updatestate!(plant, u)      # update plant simulator
    end
    res = SimResult(mpc, U_data, Y_data; Ry_data)
    return res
end

The setmodel! method must be called after solving the optimization problem with moveinput!, and before updating the state estimate with updatestate!. The SimResult object is for plotting purposes only. The adaptive LinMPC performances are similar to the nonlinear MPC, both for the 180° setpoint:

res_slin = test_slmpc(model, mpc3, [180], plant, x_0=[0, 0])
plot(res_slin)

plot10_NonLinMPC

and the 10° step disturbance:

res_slin = test_slmpc(model, mpc3, [180], plant, x_0=[π, 0], y_step=[10])
plot(res_slin)

plot11_NonLinMPC

The computations of the successive linearization MPC are about 125 times faster than the nonlinear MPC (0.00012 s per time steps versus 0.015 s per time steps, on average), an impressive gain for similar closed-loop performances!

  • 1Arnström, D., Bemporad, A., and Axehill, D. (2022). A dual active-set solver for embedded quadratic programming using recursive LDLᵀ updates. IEEE Trans. Autom. Contr., 67(8). https://doi.org/doi:10.1109/TAC.2022.3176430.