Skip to content

Commit 476552c

Browse files
authored
Merge pull request #36 from JuliaControl/h2syn
improvements to H2 design
2 parents 9909007 + b52cd7b commit 476552c

11 files changed

Lines changed: 225 additions & 77 deletions

examples/hinf_example_tank.jl

Lines changed: 7 additions & 17 deletions
Original file line numberDiff line numberDiff line change
@@ -2,14 +2,10 @@ using ControlSystems, RobustAndOptimalControl
22
using Plots
33
using LinearAlgebra
44
"""
5-
This is a simple SISO example with integrator dynamics corresponding to the
6-
quad tank process in the lab.
7-
8-
The example can be set to visualize and save plots using the variables
9-
makeplots - true/false (true if plots are to be generated, false for testing)
10-
SavePlots - true/false (true if plots are to be saved, false for testing)
5+
This is a simple SISO example with integrator dynamics corresponding to a
6+
quad tank process.
117
"""
12-
makeplots = true
8+
139

1410
# Define the proces parameters
1511
k1, k2, kc, g = 3.33, 3.35, 0.5, 981
@@ -57,16 +53,10 @@ C, γ = hinfsynthesize(P)
5753
Pcl, S, CS, T = hinfsignals(P, G, C)
5854

5955
## Plot the specifications
60-
# TODO figure out why I get segmentation errors when using ss instead of tf for
61-
# the weighting functions, makes no sense at all
62-
if makeplots
63-
specificationplot([S, CS, T], [WS[1,1], 0.01, WT[1,1]], γ)
56+
specificationplot([S, CS, T], [WS[1,1], 0.01, WT[1,1]], γ)
6457

6558
## Plot the closed loop gain from w to z
66-
# TODO figure out why the legends don't seem to work in this case
67-
68-
specificationplot(Pcl, γ; s_labels=["\$\\sigma(P_{cl}(j\\omega))\$"], w_labels=["\$\\gamma\$"])
59+
specificationplot(Pcl, γ; s_labels=["\$\\sigma(P_{cl}(j\\omega))\$"], w_labels=["\$\\gamma\$"])
6960

70-
times = [i for i in range(0, stop=300, length=10000)]
71-
plot(step(T, times))
72-
end
61+
times = [i for i in range(0, stop=300, length=10000)]
62+
plot(step(T, times))

src/ExtendedStateSpace.jl

Lines changed: 5 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -179,8 +179,13 @@ function Base.getproperty(sys::ExtendedStateSpace, s::Symbol)
179179
end
180180
end
181181

182+
Base.propertynames(sys::ExtendedStateSpace) = (:A, :B, :C, :D, :B1, :B2, :C1, :C2, :D11, :D12, :D21, :D22, :Ts, :timeevol, :nx, :ny, :nu, :nw, :nz, :zinds, :yinds, :winds, :uinds)
183+
182184
ControlSystems.StateSpace(s::ExtendedStateSpace) = ss(ssdata(s)..., s.timeevol)
183185

186+
"""
187+
A, B1, B2, C1, C2, D11, D12, D21, D22 = ssdata_e(sys)
188+
"""
184189
ssdata_e(sys::ExtendedStateSpace) = sys.A,
185190
sys.B1,
186191
sys.B2,

src/glover_mcfarlane.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -154,8 +154,8 @@ end
154154
155155
For discrete systems, the `info` tuple contains also feedback gains `F, L` and observer gain `Hkf` such that the controller on observer form is given by
156156
```math
157-
x^+ = Ax + Bu + H_{kf}*(Cx - y)\\\\
158-
u = Fx + L*(Cx - y)
157+
x^+ = Ax + Bu + H_{kf} (Cx - y)\\\\
158+
u = Fx + L (Cx - y)
159159
```
160160
Note, this controller is *not* strictly proper, i.e., it has a non-zero D matrix.
161161
The controller can be transformed to observer form for the scaled plant (`info.Gs`)

src/h2_design.jl

Lines changed: 30 additions & 15 deletions
Original file line numberDiff line numberDiff line change
@@ -10,12 +10,18 @@ If `γ = nothing`, use the formulas for H₂ in Ch 14.5. If γ is a large value,
1010
function h2synthesize(P::ExtendedStateSpace, γ = nothing)
1111

1212
if γ === nothing
13+
iszero(P.D11) && iszero(P.D22) || error("D11 and D22 must be zero for the standard H2 formulation to be used, try calling with γ=1000 to use the H∞ formulation.")
1314
X2, Y2, F2, L2 = _solvematrixequations2(P)
1415
Â2 = P.A + P.B2*F2 + L2*P.C2
15-
K = ss(Â2, -L2, F2, 0)
16+
# Af2 = P.A + P.B2*F2
17+
# C1f2 = P.C1 + P.D12*F2
18+
# Gc = ss(Af2, I(size(Af2, 1)), C1f2, 0)
19+
K = ss(Â2, -L2, F2, 0, P.timeevol)
1620
return K, lft(ss(P), K)
1721
end
1822

23+
iscontinuous(P) || throw(ArgumentError("h2syn with specified γ is only supported for continuous systems."))
24+
1925
P̄, Ltrans12, Rtrans12, Ltrans21, Rtrans21 = _transformp2pbar(P)
2026
X2, Y2, F2, L2 = _solvematrixequations(P̄, γ)
2127

@@ -52,26 +58,35 @@ function _solvematrixequations2(P::ExtendedStateSpace)
5258
B2 = P.B2
5359
C1 = P.C1
5460
C2 = P.C2
55-
D11 = P.D11
5661
D12 = P.D12
5762
D21 = P.D21
58-
D22 = P.D22
5963

60-
P1 = size(C1, 1)
61-
P2 = size(C2, 1)
62-
M1 = size(B1, 2)
63-
M2 = size(B2, 2)
64+
# P1 = size(C1, 1)
65+
# P2 = size(C2, 1)
66+
# M1 = size(B1, 2)
67+
# M2 = size(B2, 2)
68+
69+
# HX = [A zeros(size(A)); -C1'*C1 -A'] - [B2; -C1'*D12] * [D12'*C1 B2']
70+
# HY = [A' zeros(size(A)); -B1*B1' -A] - [C2'; -B1*D21'] * [D21*B1' C2]
71+
72+
# # Solve matrix equations
73+
# X2 = _solvehamiltonianare(HX)
74+
# Y2 = _solvehamiltonianare(HY)
75+
76+
# # These formulas appear under ch14.5, but different formulas where D12'C1=0 and B1 * D21'=0 appear
77+
# # in ch 14.9.1, the difference is explained in ch13. See also LQGProblem(P::ExtendedStateSpace)
78+
# F2 = (B2'X2 + D12'C1)
79+
# L2 = (B1 * D21' + Y2 * C2')
6480

65-
HX = [A zeros(size(A)); -C1'*C1 -A'] - [B2; -C1'*D12] * [D12'*C1 B2']
66-
HY = [A' zeros(size(A)); -B1*B1' -A] - [C2'; -B1*D21'] * [D21*B1' C2]
6781

68-
# Solve matrix equations
69-
X2 = _solvehamiltonianare(HX)
70-
Y2 = _solvehamiltonianare(HY)
82+
# This is more robust than the approach above
83+
# https://ocw.snu.ac.kr/sites/default/files/NOTE/3938.pdf eq 73-74
7184

72-
F2 = -(B2'X2 + D12'C1)
85+
fun = iscontinuous(P) ? MatrixEquations.arec : MatrixEquations.ared
7386

74-
L2 = -(B1 * D21' + Y2 * C2')
87+
X2, _, F2 = fun(A, B2, D12'D12, C1'C1, C1'D12)
88+
Y2, _, L2 = fun(A', C2', D21*D21', B1*B1', B1*D21')
89+
L2 = L2'
7590

76-
return X2, Y2, F2, L2
91+
return X2, Y2, -F2, -L2
7792
end

src/hinfinity_design.jl

Lines changed: 35 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -254,15 +254,16 @@ function _synthesizecontroller(
254254
D1122 = D11[(P1-M2+1):P1, (M1-P2+1):M1]
255255

256256
# Equation 19
257-
D11hat = ((-D1121 * D1111') / (γ² * I - D1111 * D1111')) * D1112 - D1122
257+
J1 = (γ² * I - D1111 * D1111')
258+
D11hat = ((-D1121 * D1111') / J1) * D1112 - D1122
258259

259260
# Equation 20
260261
D12hatD12hat = I - (D1121 / (γ² * I - D1111' * D1111)) * D1121'
261262
_assertrealandpsd(D12hatD12hat; msg = " in equation (20)")
262263
D12hat = cholesky(Hermitian(D12hatD12hat)).L
263264

264265
# Equation 21
265-
D21hatD21hat = I - (D1112' / (γ² * I - D1111 * D1111')) * D1112
266+
D21hatD21hat = I - (D1112' / J1) * D1112
266267
_assertrealandpsd(D21hatD21hat; msg = " in equation (21)")
267268
D21hat = cholesky(Hermitian(D21hatD21hat)).U
268269

@@ -306,6 +307,29 @@ function _synthesizecontroller(
306307
return ss(Ac, Bc[:, 1:P2], Cc[1:M2, :], D11c)
307308
end
308309

310+
"""
311+
rqr(D, γ=1)
312+
313+
"Regularized" qr factorization. This struct represents \$(D'D + γI)\$ without forming \$D'D\$
314+
Note: this does not support negative γ or \$(γI - D'D)\$
315+
Supported operations: `\\,/,*`, i.e., it behaves also like a matrix (unlike the standard `QR` factorization object).
316+
"""
317+
struct rqr{T, PT, DGT}
318+
P::PT
319+
DG::DGT
320+
function rqr(D, γ=1)
321+
P = (D'D) + γ*I
322+
DG = qr([D; (γ)I])
323+
new{eltype(P), typeof(P), typeof(DG)}(P, DG)
324+
end
325+
end
326+
327+
Base.:\(d::rqr, b) = (d.DG.R\(adjoint(d.DG.R)\b))
328+
Base.:/(b, d::rqr) = ((b/d.DG.R)/adjoint(d.DG.R))
329+
Base.:*(d::rqr, b) = (d.P*b)
330+
Base.:*(b, d::rqr) = (b*d.P)
331+
332+
309333
"""
310334
_assertrealandpsd(A::AbstractMatrix, msg::AbstractString)
311335
@@ -374,7 +398,7 @@ function _solvehamiltonianare(H)#::AbstractMatrix{<:LinearAlgebra.BlasFloat})
374398
U11 = S.Z[1:div(m, 2), 1:div(n, 2)]
375399
U21 = S.Z[div(m, 2)+1:m, 1:div(n, 2)]
376400

377-
return U21 / U11
401+
return U21 * pinv(U11)
378402
end
379403

380404
"""
@@ -447,11 +471,12 @@ function _γiterations(
447471
)
448472

449473
T = typeof(P.A)
474+
ET = eltype(T)
450475
XinfFeasible, YinfFeasible, FinfFeasible, HinfFeasible, gammFeasible =
451476
T(undef,0,0), T(undef,0,0), T(undef,0,0), T(undef,0,0), nothing
452477

453-
gl, gu = interval
454-
gl = max(1e-3, gl)
478+
gl, gu = ET.(interval)
479+
gl = max(ET(1e-3), gl)
455480
iters = ceil(Int, log2((gu-gl+1e-16)/gtol))
456481

457482
for iteration = 1:iters
@@ -594,6 +619,10 @@ can be both MIMO and SISO, both in tf and ss forms). Valid inputs for the
594619
weighting functions are empty arrays, numbers (static gains), and `LTISystem`s.
595620
"""
596621
function hinfpartition(G, WS, WU, WT)
622+
WS isa LTISystem && common_timeevol(G,WS)
623+
WU isa LTISystem && common_timeevol(G,WU)
624+
WT isa LTISystem && common_timeevol(G,WT)
625+
te = G.timeevol
597626
# # Convert the systems into state-space form
598627
Ag, Bg, Cg, Dg = _input2ss(G)
599628
Asw, Bsw, Csw, Dsw = _input2ss(WS)
@@ -726,7 +755,7 @@ function hinfpartition(G, WS, WU, WT)
726755
Dyw = Matrix{Float64}(I, mCg, nDuw)
727756
Dyu = -Dg
728757

729-
P = ss(A, Bw, Bu, Cz, Cy, Dzw, Dzu, Dyw, Dyu)
758+
P = ss(A, Bw, Bu, Cz, Cy, Dzw, Dzu, Dyw, Dyu, te)
730759

731760
end
732761

0 commit comments

Comments
 (0)