% Optimal Control % ๐ค Sรฉbastien Boisgรฉrault
-
๐ Documents (GitHub)
-
ยฉ๏ธ License CC BY 4.0
๐ | Code | ๐ | Worked Example |
๐ | Graph | ๐งฉ | Exercise |
๐ท๏ธ | Definition | ๐ป | Numerical Method |
๐ | Theorem | ๐งฎ | Analytical Method |
๐ | Remark | ๐ง | Theory |
โน๏ธ | Information | ๐๏ธ | Hint |
Warning | ๐ | Solution |
from numpy import *
from numpy.linalg import *
from matplotlib.pyplot import *
from scipy.integrate import solve_ivp
from scipy.linalg import solve_continuous_are
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
# Python 3.x Standard Library
import gc
import os
# Third-Party Packages
import numpy as np; np.seterr(all="ignore")
import numpy.linalg as la
import scipy.misc
import matplotlib as mpl; mpl.use("Agg")
import matplotlib.pyplot as pp
import matplotlib.axes as ax
import matplotlib.patches as pa
#
# Matplotlib Configuration & Helper Functions
# --------------------------------------------------------------------------
# TODO: also reconsider line width and markersize stuff "for the web
# settings".
fontsize = 10
width = 345 / 72.27
height = width / (16/9)
rc = {
"text.usetex": True,
"pgf.preamble": r"\usepackage{amsmath,amsfonts,amssymb}",
#"font.family": "serif",
"font.serif": [],
#"font.sans-serif": [],
"legend.fontsize": fontsize,
"axes.titlesize": fontsize,
"axes.labelsize": fontsize,
"xtick.labelsize": fontsize,
"ytick.labelsize": fontsize,
"figure.max_open_warning": 100,
#"savefig.dpi": 300,
#"figure.dpi": 300,
"figure.figsize": [width, height],
"lines.linewidth": 1.0,
}
mpl.rcParams.update(rc)
# Web target: 160 / 9 inches (that's ~45 cm, this is huge) at 90 dpi
# (the "standard" dpi for Web computations) gives 1600 px.
width_in = 160 / 9
def save(name):
cwd = os.getcwd()
root = os.path.dirname(os.path.realpath(__file__))
os.chdir(root)
pp.savefig(name + ".svg")
os.chdir(cwd)
def set_ratio(ratio=1.0, bottom=0.1, top=0.1, left=0.1, right=0.1):
height_in = (1.0 - left - right)/(1.0 - bottom - top) * width_in / ratio
pp.gcf().set_size_inches((width_in, height_in))
pp.gcf().subplots_adjust(bottom=bottom, top=1.0-top, left=left, right=1.0-right)
width
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
Limitations of Pole Assignment
-
It is not always obvious what set of poles we should target (especially for large systems),
-
We do not control explicitly the trade-off between "speed of convergence" and "intensity of the control" (large input values maybe costly or impossible).
Let
where
-
$A \in \mathbb{R}^{n\times n}$ ,$B \in \mathbb{R}^{m\times n}$ and -
$x(0) = x_0 \in \mathbb{R}^n$ is given.
Find
where:
-
$Q \in \mathbb{R}^{n \times n}$ and$R \in \mathbb{R}^{m\times m}$ , -
(to be continued ...)
-
$Q$ and$R$ are symmetric ($R^t = R$ and$Q^t = Q$ ), -
$Q$ and$R$ are positive definite (denoted "$>0$")$$x^t Q x \geq 0 , \mbox{ and } , x^t Q x = 0 , \mbox{ iff }, x=0$$ and
$$u^t R u \geq 0 , \mbox{ and } , u^t R u = 0 , \mbox{ iff }, u=0.$$
If
with
When we minimize
-
Only the relative values of
$q$ and$r$ matters. -
Large values of
$q$ penalize strongly non-zero states:$\Rightarrow$ fast convergence. -
Large values of
$r$ penalize strongly non-zero inputs:$\Rightarrow$ small input values.
If
with
Thus we can control the cost of each component of
Assume that
-
There is an optimal solution; it is a linear feedback
$$u = - K x$$ -
The closed-loop dynamics is asymptotically stable.
-
The gain matrix
$K$ is given by$$ K = R^{-1} B^t \Pi, $$
where
$\Pi \in \mathbb{R}^{n \times n}$ is the unique matrix such that$\Pi^t = \Pi$ ,$\Pi > 0$ and$$ \Pi B R^{-1} B^t \Pi - \Pi A - A^t \Pi - Q = 0. $$
Consider the double integrator
\left[\begin{array}{cx} 0 & 1 \ 0 & 0\end{array}\right] \left[\begin{array}{c} x \ \dot{x} \end{array}\right] + \left[\begin{array}{c} 0 \ 1 \end{array}\right] u $$
(in standard form)
A = array([[0, 1], [0, 0]])
B = array([[0], [1]])
Q = array([[1, 0], [0, 1]])
R = array([[1]])
Pi = solve_continuous_are(A, B, Q, R)
K = inv(R) @ B.T @ Pi
It is stable:
eigenvalues, _ = eig(A - B @ K)
assert all([real(s) < 0 for s in eigenvalues])
figure()
x = [real(s) for s in eigenvalues]
y = [imag(s) for s in eigenvalues]
plot(x, y, "kx")
grid(True)
title("Eigenvalues")
axis("square")
axis([-5, 5, -5, 5])
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
save("images/poles-LQ")
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::: slides :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
y0 = [1.0, 1.0]
def f(t, x):
return (A - B @ K) @ x
result = solve_ivp(
f, t_span=[0, 10], y0=y0, max_step=0.1
)
t = result["t"]
x1 = result["y"][0]
x2 = result["y"][1]
u = - (K @ result["y"]).flatten() # vect. -> scalar
figure()
plot(t, x1, label="$x_1$")
plot(t, x2, label="$x_2$")
plot(t, u, label="$u$")
xlabel("$t$")
grid(True)
legend(loc="lower right")
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
tight_layout()
save("images/poles-LQ-traj")
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::: slides :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
Q = array([[10, 0], [0, 10]])
R = array([[1]])
Pi = solve_continuous_are(A, B, Q, R)
K = inv(R) @ B.T @ Pi
eigenvalues, _ = eig(A - B @ K)
assert all([real(s) < 0 for s in eigenvalues])
figure()
x = [real(s) for s in eigenvalues]
y = [imag(s) for s in eigenvalues]
plot(x, y, "kx")
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
xticks(arange(-5, 6)); yticks(arange(-5, 6))
plot([0, 0], [-5, 5], "k")
plot([-5, 5], [0, 0], "k")
grid(True)
title("Eigenvalues")
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
axis("square")
axis([-5, 5, -5, 5])
save("images/poles-LQ-2")
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::: slides :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
result = solve_ivp(
f, t_span=[0, 10], y0=y0, max_step=0.1
)
t = result["t"]
x1 = result["y"][0]
x2 = result["y"][1]
u = - (K @ result["y"]).flatten() # vect. -> scalar
figure()
plot(t, x1, label="$x_1$")
plot(t, x2, label="$x_2$")
plot(t, u, label="$u$")
xlabel("$t$")
grid(True)
legend(loc="lower right")
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
tight_layout()
save("images/poles-LQ-2-traj")
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::: slides :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
Q = array([[1, 0], [0, 1]])
R = array([[10]])
Pi = solve_continuous_are(A, B, Q, R)
K = inv(R) @ B.T @ Pi
eigenvalues, _ = eig(A - B @ K)
assert all([real(s) < 0 for s in eigenvalues])
figure()
x = [real(s) for s in eigenvalues]
y = [imag(s) for s in eigenvalues]
plot(x, y, "kx")
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
xticks(arange(-5, 6)); yticks(arange(-5, 6))
plot([0, 0], [-5, 5], "k")
plot([-5, 5], [0, 0], "k")
grid(True)
title("Eigenvalues")
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
axis("square")
axis([-5, 5, -5, 5])
save("images/poles-LQ-3")
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::: slides :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
result = solve_ivp(
f, t_span=[0, 10], y0=y0, max_step=0.1
)
t = result["t"]
x1 = result["y"][0]
x2 = result["y"][1]
u = - (K @ result["y"]).flatten() # vect. -> scalar
figure()
plot(t, x1, label="$x_1$")
plot(t, x2, label="$x_2$")
plot(t, u, label="$u$")
xlabel("$t$")
grid(True)
legend(loc="lower right")
::: hidden :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
tight_layout()
save("images/poles-LQ-3-traj")
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::: slides :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
Consider the controllable dynamics
and
Let
Show that
What is the value of
We know that
Since
and thus
Since
Since the system is controllable, the optimal control makes the
origin of the closed-loop system asymptotically stable. Consequently,