Skip to content

Commit 3b1a46d

Browse files
authored
add root-locus methods for matrix gains (#1012)
* add root-locus methods for matrix gains * rm comment * add kwargs...
1 parent 948aa79 commit 3b1a46d

File tree

3 files changed

+224
-6
lines changed

3 files changed

+224
-6
lines changed

src/root_locus.jl

Lines changed: 121 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -3,7 +3,7 @@ import ControlSystemsBase.RootLocusResult
33
@userplot Rlocusplot
44

55

6-
function getpoles(G, K::Number)
6+
function getpoles(G, K::Number; kwargs...)
77
issiso(G) || error("root locus only supports SISO systems")
88
G isa TransferFunction || (G = tf(G))
99
P = numpoly(G)[]
@@ -69,17 +69,115 @@ function getpoles(G, K::AbstractVector{T}) where {T<:Number}
6969
copy(poleout'), K
7070
end
7171

72+
"""
73+
getpoles(sys::StateSpace, K::AbstractMatrix; tol = 1e-2, initial_stepsize = 1e-3, output=false)
74+
75+
Compute the poles of the closed-loop system defined by `sys` with feedback gains `γ*K` where `γ` is a scalar that ranges from 0 to 1.
76+
77+
If `output = true`, `K` is assumed to be an output feedback matrix of dim `(nu, ny)`
78+
"""
79+
function getpoles(sys::StateSpace, K_matrix::AbstractMatrix; tol = 1e-2, initial_stepsize = 1e-3, output=false)
80+
(; A, B, C) = sys
81+
nx = size(A, 1) # State dimension
82+
ny = size(C, 1) # Output dimension
83+
tol = tol*nx # Scale tolerance with state dimension
84+
# Check for compatibility of K_matrix dimensions with B
85+
if size(K_matrix, 2) != (output ? ny : nx)
86+
error("The number of columns in K_matrix ($(size(K_matrix, 2))) must match the state dimension ($(nx)) or output dimension ($(ny)) depending on whether output feedback is used.")
87+
end
88+
if size(K_matrix, 1) != size(B, 2)
89+
error("The number of rows in K_matrix ($(size(K_matrix, 1))) must match the number of inputs (columns of B, which is $(size(B, 2))).")
90+
end
91+
92+
if output
93+
# We bake C into K here to avoid repeated multiplications below
94+
K_matrix = K_matrix * C
95+
end
96+
97+
poleout_list = Vector{Vector{ComplexF64}}() # To store pole sets at each accepted step
98+
k_scalars_collected = Float64[] # To store accepted k_scalar values
99+
100+
prevpoles = ComplexF64[] # Initialize prevpoles for the first iteration
101+
102+
stepsize = initial_stepsize
103+
k_scalar = 0.0
104+
105+
# Initial poles at k_scalar = 0.0
106+
A_cl_initial = A - 0.0 * B * K_matrix
107+
initial_poles = eigvals(A_cl_initial)
108+
push!(poleout_list, initial_poles)
109+
push!(k_scalars_collected, 0.0)
110+
prevpoles = initial_poles # Set prevpoles for the first actual step
111+
112+
while k_scalar < 1.0
113+
# Propose a new k_scalar value
114+
next_k_scalar = min(1.0, k_scalar + stepsize)
115+
116+
# Calculate poles for the proposed next_k_scalar
117+
A_cl_proposed = A - next_k_scalar * B * K_matrix
118+
current_poles_proposed = eigvals(A_cl_proposed)
119+
120+
# Calculate cost using Hungarian algorithm
121+
D = zeros(nx, nx)
122+
for r in 1:nx
123+
for c in 1:nx
124+
D[r, c] = abs(current_poles_proposed[r] - prevpoles[c])
125+
end
126+
end
127+
assignment, cost = Hungarian.hungarian(D)
128+
129+
# Adaptive step size logic
130+
if cost > 2 * tol # Cost is too high, reject step and reduce stepsize
131+
stepsize /= 2.0
132+
# Ensure stepsize doesn't become too small
133+
if stepsize < 100eps()
134+
@warn "Step size became extremely small, potentially stuck. Breaking loop."
135+
break
136+
end
137+
# Do not update k_scalar, try again with smaller stepsize
138+
else # Step is acceptable or too small
139+
# Sort poles using the assignment from Hungarian algorithm
140+
temppoles = zeros(ComplexF64, nx)
141+
for j = 1:nx
142+
temppoles[assignment[j]] = current_poles_proposed[j]
143+
end
144+
current_poles_sorted = temppoles
145+
146+
# Accept the step
147+
push!(poleout_list, current_poles_sorted)
148+
push!(k_scalars_collected, next_k_scalar)
149+
prevpoles = current_poles_sorted # Update prevpoles for the next iteration
150+
k_scalar = next_k_scalar # Advance k_scalar
151+
152+
if cost < tol # Cost is too low, increase stepsize
153+
stepsize *= 1.1
154+
end
155+
# Cap stepsize to prevent overshooting 1.0 significantly in a single step
156+
stepsize = min(stepsize, 1e-1)
157+
end
158+
159+
# Break if k_scalar has reached or exceeded 1.0
160+
if k_scalar >= 1.0
161+
break
162+
end
163+
end
164+
165+
return hcat(poleout_list...)' |> copy, k_scalars_collected .* Ref(K_matrix) # Return transposed pole matrix and k_values
166+
end
167+
72168

73169
"""
74170
roots, Z, K = rlocus(P::LTISystem, K = 500)
75171
76172
Compute the root locus of the SISO LTISystem `P` with a negative feedback loop and feedback gains between 0 and `K`. `rlocus` will use an adaptive step-size algorithm to determine the values of the feedback gains used to generate the plot.
77173
78174
`roots` is a complex matrix containing the poles trajectories of the closed-loop `1+k⋅G(s)` as a function of `k`, `Z` contains the zeros of the open-loop system `G(s)` and `K` the values of the feedback gain.
175+
176+
If `K` is a matrix and `P` a `StateSpace` system, the poles are computed as `K` ranges from `0*K` to `1*K`. In this case, `K` is assumed to be a state-feedback matrix of dimension `(nu, nx)`. To compute the poles for output feedback, use, pass `output = true` and `K` of dimension `(nu, ny)`.
79177
"""
80-
function rlocus(P, K)
178+
function rlocus(P, K; kwargs...)
81179
Z = tzeros(P)
82-
roots, K = getpoles(P,K)
180+
roots, K = getpoles(P, K; kwargs...)
83181
ControlSystemsBase.RootLocusResult(roots, Z, K, P)
84182
end
85183

@@ -89,6 +187,7 @@ rlocus(P; K=500) = rlocus(P, K)
89187
# This will be called on plot(rlocus(sys, args...))
90188
@recipe function rootlocusresultplot(r::RootLocusResult)
91189
roots, Z, K = r
190+
array_K = eltype(K) <: AbstractArray
92191
redata = real.(roots)
93192
imdata = imag.(roots)
94193
all_redata = [vec(redata); real.(Z)]
@@ -103,7 +202,9 @@ rlocus(P; K=500) = rlocus(P, K)
103202
form(k, p) = Printf.@sprintf("%.4f", k) * " pole=" * Printf.@sprintf("%.3f%+.3fim", real(p), imag(p))
104203
@series begin
105204
legend --> false
106-
hover := "K=" .* form.(K,roots)
205+
if !array_K
206+
hover := "K=" .* form.(K,roots)
207+
end
107208
label := ""
108209
redata, imdata
109210
end
@@ -121,6 +222,15 @@ rlocus(P; K=500) = rlocus(P, K)
121222
label --> "Open-loop poles"
122223
redata[1,:], imdata[1,:]
123224
end
225+
if array_K
226+
@series begin
227+
seriestype := :scatter
228+
markershape --> :diamond
229+
markersize --> 10
230+
label --> "Closed-loop poles"
231+
redata[end,:], imdata[end,:]
232+
end
233+
end
124234
end
125235

126236

@@ -129,5 +239,10 @@ end
129239
130240
Plot the root locus of the SISO LTISystem `P` as computed by `rlocus`.
131241
"""
132-
@recipe rlocusplot(::Type{Rlocusplot}, p::Rlocusplot; K=500) =
133-
rlocus(p.args[1]; K=K)
242+
@recipe function rlocusplot(::Type{Rlocusplot}, p::Rlocusplot; K=500, output=false)
243+
if length(p.args) >= 2
244+
rlocus(p.args[1], p.args[2]; output)
245+
else
246+
rlocus(p.args[1]; K=K)
247+
end
248+
end

test/runtests_controlsystems.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -20,5 +20,6 @@ using ControlSystems
2020
@testset "rootlocus" begin
2121
@info "Testing rootlocus"
2222
include("test_rootlocus.jl")
23+
include("test_root_locus_matrix.jl")
2324
end
2425
end

test/test_root_locus_matrix.jl

Lines changed: 102 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,102 @@
1+
using ControlSystems
2+
using ControlSystems: getpoles
3+
using Test
4+
using Plots
5+
6+
@testset "Root Locus with Matrix Feedback" begin
7+
8+
@testset "State Feedback Matrix" begin
9+
# Test with double mass model
10+
P = DemoSystems.double_mass_model()
11+
12+
# Design a simple state feedback matrix
13+
K_state = [1.0 0.5 0.2 0.1] # nu x nx matrix (1x4)
14+
15+
# Test that getpoles works with matrix feedback
16+
roots, K_values = getpoles(P, K_state)
17+
18+
# Basic sanity checks
19+
@test size(roots, 2) == size(P.A, 1) # Should have nx poles
20+
@test length(K_values) == size(roots, 1) # K_values should match number of steps
21+
@test eltype(K_values) <: AbstractMatrix # K_values should be matrices
22+
23+
# Test that initial K_value is zero matrix
24+
@test all(K_values[1] .≈ 0.0)
25+
26+
# Test that final K_value is the input matrix
27+
@test K_values[end] K_state
28+
29+
# Test poles are continuous (no sudden jumps)
30+
# The poles should change smoothly
31+
pole_diffs = diff(roots, dims=1)
32+
max_diff = maximum(abs, pole_diffs)
33+
@test max_diff < 1.0 # Poles shouldn't jump too much between steps
34+
35+
# Test that rlocus works with matrix K
36+
result = rlocus(P, K_state)
37+
@test result isa ControlSystemsBase.RootLocusResult
38+
@test size(result.roots) == size(roots)
39+
40+
# Test plotting (should not error)
41+
plt = plot(result)
42+
@test plt isa Plots.Plot
43+
end
44+
45+
@testset "Output Feedback Matrix" begin
46+
P = DemoSystems.double_mass_model(outputs=1:2)
47+
48+
# Design output feedback matrix
49+
# P has 2 outputs, 1 input, so K should be 1x2
50+
K_output = [1.0 0.5] # nu x ny matrix (1x2)
51+
52+
# Test output feedback
53+
roots, K_values = getpoles(P, K_output; output=true)
54+
55+
# Basic checks
56+
@test size(roots, 2) == size(P.A, 1) # Should have nx poles
57+
@test length(K_values) == size(roots, 1)
58+
@test eltype(K_values) <: AbstractMatrix
59+
60+
# Test rlocus with output feedback
61+
result = rlocus(P, K_output; output=true)
62+
@test result isa ControlSystemsBase.RootLocusResult
63+
64+
# Test plotting
65+
plt = plot(result)
66+
@test plt isa Plots.Plot
67+
end
68+
69+
@testset "Error Handling" begin
70+
P = DemoSystems.double_mass_model()
71+
72+
# Test wrong dimensions for state feedback
73+
K_wrong = [1.0 0.5 0.2] # Wrong size (1x3 instead of 1x4)
74+
@test_throws ErrorException getpoles(P, K_wrong)
75+
76+
# Test wrong dimensions for output feedback
77+
K_wrong_output = [1.0 0.5 0.2] # Wrong size (1x3 instead of 1x2)
78+
@test_throws ErrorException getpoles(P, K_wrong_output; output=true)
79+
80+
# Test wrong number of rows
81+
K_wrong_rows = [1.0 0.5 0.2 0.1; 2.0 1.0 0.4 0.2] # 2x4 instead of 1x4
82+
@test_throws ErrorException getpoles(P, K_wrong_rows)
83+
end
84+
85+
86+
@testset "Adaptive Step Size" begin
87+
P = DemoSystems.double_mass_model()
88+
K_state = [1.0 0.5 0.2 0.1]
89+
90+
# Test with different tolerances
91+
roots_loose, K_loose = getpoles(P, K_state; tol=1e-1)
92+
roots_tight, K_tight = getpoles(P, K_state; tol=1e-4)
93+
94+
# Tighter tolerance should result in more steps
95+
@test length(K_tight) >= length(K_loose)
96+
97+
# Both should start and end at the same points
98+
@test roots_loose[1, :] roots_tight[1, :] # Same initial poles
99+
@test roots_loose[end, :] roots_tight[end, :] atol=1e-3 # Similar final poles
100+
end
101+
102+
end

0 commit comments

Comments
 (0)