RobustAndOptimalControl.jl

Installation

pkg> add RobustAndOptimalControl

Named systems

named_ss

Connecting systems together

connect

The following complicated feedback interconnection

                 yF
              ┌────────────────────────────────┐
              │                                │
    ┌───────┐ │  ┌───────┐ yR   ┌─────────┐    │    ┌───────┐
uF  │       │ │  │       ├──────►         │ yC │  uP│       │    yP
────►   F   ├─┴──►   R   │      │    C    ├────+────►   P   ├────┬────►
    │       │    │       │   ┌──►         │         │       │    │
    └───────┘    └───────┘   │  └─────────┘         └───────┘    │
                             │                                   │
                             └───────────────────────────────────┘

can be created by

F = named_ss(ssrand(1, 1, 2, proper=true), x=:xF, u=:uF, y=:yF)
R = named_ss(ssrand(1, 1, 2, proper=true), x=:xR, u=:uR, y=:yR)
C = named_ss(ssrand(1, 1, 2, proper=true), x=:xC, u=:uC, y=:yC)
P = named_ss(ssrand(1, 1, 3, proper=true), x=:xP, u=:uP, y=:yP)

addP = sumblock("uP = yF + yC") # Sum node before P
addC = sumblock("uC = yR - yP") # Sum node before C

connections = [
    :yP => :yP # Output to input
    :uP => :uP
    :yC => :yC
    :yF => :yF
    :yF => :uR
    :uC => :uC
    :yR => :yR
]
w1 = [:uF] # External inputs

G = connect([F, R, C, P, addP, addC], connections; w1)

If an external input is to be connected to multiple points, use a splitter to split up the signal into a set of unique names which are then used in the connections.

Uncertainty modeling

TODO: add examples

Model augmentation

TODO

$H_\infty$ and $H_2$ design

TODO Examples are available in the example folder.

LQG design

TODO

Structured singular value and diskmargin

TODO

Closed-loop analysis

Exported functions and types

Index

RobustAndOptimalControl.LQGProblemType
G = LQG(sys::AbstractStateSpace, Q1, Q2, R1, R2; qQ=0, qR=0, M = I, N = I)

Return an LQG object that describes the closed control loop around the process sys=ss(A,B,C,D) where the controller is of LQG-type. The controller is specified by weight matrices Q1,Q2 that penalizes state deviations and control signal variance respectively, and covariance matrices R1,R2 which specify state drift and measurement covariance respectively.

qQ and qR can be set to incorporate loop transfer recovery, i.e.,

L = lqr(A, B, Q1+qQ*C'C, Q2)
K = kalman(A, C, R1+qR*B*B', R2)

Increasing qQ will add more cost in output direction, e.g., encouraging the use of cheap control, while increasing qR adds fictious dynamics noise, makes the observer faster.

M is a matrix that defines the controlled variables z, i.e., the variables for which you provide reference signals. If no M is provided, the default is to consider all state variables of the system as controlled. The definitions of z and y are given below

y = C*x
z = M*x

size(M, 1) determines the size of the Q1 matrix you need to supply.

N is a matrix that defines how the dynamics noise v enters the system, i.e. If no N is provided, the default is to consider all state variables being affected by independent noise components. The definition of v is given below

x′ = A*x + B*u + N*v

size(N, 2) determines the size of the R1 matrix you need to supply.

Example

s = tf("s")
P = [1/(s+1) 2/(s+2); 1/(s+3) 1/(s-1)]
sys = ss(P)
eye(n) = Matrix{Float64}(I,n,n) # For convinience

qQ = 1
qR = 1
Q1 = 10eye(4)
Q2 = 1eye(2)
R1 = 1eye(6)
R2 = 1eye(2)

G = LQG(sys, Q1, Q2, R1, R2, qQ=qQ, qR=qR)

Gcl = G.cl
T = G.T
S = G.S
sigmaplot([S,T],exp10.(range(-3, stop=3, length=1000)))
stepplot(Gcl)
RobustAndOptimalControl.LQGProblemMethod
LQGProblem(P::ExtendedStateSpace)

If only an ExtendedStateSpace system is provided, the system P is assumed to correspond to the H₂ optimal control problem with

C1'C1    = Q1
D12'D12  = Q2

B1*B1'   = R1
D21*D21' = R2

and an LQGProblem with the above covariance matrices is returned. The system description in the returned LQGProblem will have B1 = C1 = I. See Ch. 13 in Robust and optimal control for reference.

RobustAndOptimalControl.add_disturbanceMethod
add_disturbance(sys::AbstractStateSpace, Ad::AbstractMatrix, Cd::AbstractMatrix)

See CCS pp. 144

Arguments:

  • sys: System to augment
  • Ad: The dynamics of the disturbance
  • Cd: How the disturbance states affect the states of sys. This matrix has the shape (sys.nx, size(Ad, 1))

See also add_low_frequency_disturbance, add_resonant_disturbance

RobustAndOptimalControl.add_input_differentiatorFunction
add_input_differentiator(sys::AbstractStateSpace, ui = 1:sys.nu; goodwin=false)

Augment the output of sys with the difference u(k+1)-u(k)

Arguments:

  • ui: An index or vector of indices indicating which inputs to differentiate.
  • goodwin: If true, the difference operator will use the Goodwin δ operator, i.e., (u(k+1)-u(k)) / sys.Ts.

The augmented system will have the matrices

[A 0; 0 0]  [B; I]  [C 0; 0 -I]  [D; I]

with length(ui) added states and outputs.

RobustAndOptimalControl.add_output_differentiatorFunction
add_differentiator(sys::AbstractStateSpace{<:Discrete})

Augment the output of sys with the numerical difference (discrete-time derivative) of output, i.e., y_aug = [y; (y-y_prev)/sys.Ts] To add both an integrator and a differentiator to a SISO system, use

RobustAndOptimalControl.add_output_integratorFunction

addoutputintegrator(sys::AbstractStateSpace{<:Discrete}, ind = 1; ϵ = 0)

Augment the output of sys with the integral of output at index ind, i.e., y_aug = [y; ∫y[ind]] To add both an integrator and a differentiator to a SISO system, use

Gd = add_output_integrator(add_output_differentiator(G), 1)

Note: numerical integration is subject to numerical drift. If the output of the system corresponds to, e.g., a velocity reference and the integral to position reference, consider methods for mitigating this drift.

RobustAndOptimalControl.bilinearc2dMethod
bilinearc2d(Ac::AbstractArray, Bc::AbstractArray, Cc::AbstractArray, Dc::AbstractArray, Ts::Number; tolerance=1e-12)

Balanced Bilinear transformation in State-Space. This method computes a discrete time equivalent of a continuous-time system, such that

G_d(z) = s2z[G_c(s)]

in a manner which accomplishes the following (i) Preserves the infinity L-infinity norm over the transformation (ii) Finds a system which balances B and C, in the sense that ||B||2=||C||2 (iii) Satisfies Gc(s) = z2s[s2z[Gc(s)]] for some map z2s[]

RobustAndOptimalControl.bilinearc2dMethod
bilinearc2d(sys::ExtendedStateSpace, Ts::Number)

Applies a Balanced Bilinear transformation to a discrete-time extended statespace object

RobustAndOptimalControl.bilineard2cMethod
bilineard2c(Ad::AbstractArray, Bd::AbstractArray, Cd::AbstractArray, Dd::AbstractArray, Ts::Number; tolerance=1e-12)

Balanced Bilinear transformation in State-Space. This method computes a continuous time equivalent of a discrete time system, such that

G_c(z) = z2s[G_d(z)]

in a manner which accomplishes the following (i) Preserves the infinity L-infinity norm over the transformation (ii) Finds a system which balances B and C, in the sense that ||B||2=||C||2 (iii) Satisfies Gd(z) = s2z[z2s[Gd(z)]] for some map s2z[]

RobustAndOptimalControl.closedloopFunction
closedloop(l::LQGProblem, L = lqr(l), K = kalman(l))

Closed-loop system as defined in Glad and Ljung eq. 8.28

The return value will be the closed loop from reference only, other disturbance signals (B1) are ignored. See feedback for a more advanced option.

RobustAndOptimalControl.comp_sensitivityMethod

See output_comp_sensitivity

         ▲
         │e₁
         │  ┌─────┐
d₁────+──┴──►  P  ├─────┬──►e₄
      │-    └─────┘     │
      │                 │
      │     ┌─────┐     │
 e₂◄──┴─────┤  C  ◄──┬──+───d₂
            └─────┘  │
                     │e₃
                     ▼
  • Input sensitivity is the transfer function from d₁ to e₁, (I + CP)⁻¹
  • Output sensitivity is the transfer function from d₂ to e₂, (I + PC)⁻¹
  • Input complementary sensitivity is the transfer function from d₁ to e₂, CP(I + CP)⁻¹
  • Output complementary sensitivity is the transfer function from d₂ to e₄, PC(I + PC)⁻¹
RobustAndOptimalControl.connectMethod
connect(systems, connections; w1, z1 = (:), verbose = true, kwargs...)

Create complicated feedback interconnection.

Addition and subtraction nodes are achieved by creating a linear combination node, i.e., a system with a D matrix only.

Arguments:

  • systems: A vector of named systems to be connected
  • connections: a vector of pairs indicating output => input mappings.
    • u1: input mappings (alternative input argument)
    • y1: output mappings (alternative input argument)
  • w1: external signals
  • z1: outputs (can overlap with y1)
  • verbose: Issue warnings for signals that have no connection

Example: The following complicated feedback interconnection

                 yF
              ┌────────────────────────────────┐
              │                                │
    ┌───────┐ │  ┌───────┐ yR   ┌─────────┐    │    ┌───────┐
uF  │       │ │  │       ├──────►         │ yC │  uP│       │    yP
────►   F   ├─┴──►   R   │      │    C    ├────+────►   P   ├────┬────►
    │       │    │       │   ┌──►         │         │       │    │
    └───────┘    └───────┘   │  └─────────┘         └───────┘    │
                             │                                   │
                             └───────────────────────────────────┘

can be created by

F = named_ss(ssrand(1, 1, 2, proper=true), x=:xF, u=:uF, y=:yF)
R = named_ss(ssrand(1, 1, 2, proper=true), x=:xR, u=:uR, y=:yR)
C = named_ss(ssrand(1, 1, 2, proper=true), x=:xC, u=:uC, y=:yC)
P = named_ss(ssrand(1, 1, 3, proper=true), x=:xP, u=:uP, y=:yP)

addP = sumblock("uP = yF + yC") # Sum node before P
addC = sumblock("uC = yR - yP") # Sum node before C

connections = [
    :yP => :yP # Output to input
    :uP => :uP
    :yC => :yC
    :yF => :yF
    :yF => :uR
    :uC => :uC
    :yR => :yR
]
w1 = [:uF] # External inputs

G = connect([F, R, C, P, addP, addC], connections; w1)

If an external input is to be connected to multiple points, use a splitter to split up the signal into a set of unique names which are then used in the connections.

RobustAndOptimalControl.controller_reductionFunction
controller_reduction(P, K, r, out=false)

Minimize ||(K-Kᵣ) W||∞ if out=false ||W (K-Kᵣ)||∞ if out=true See Robust and Optimal Control Ch 19.1 out indicates if the weight will be applied as output or input weight.

RobustAndOptimalControl.expand_symbolMethod
expand_symbol(s::Symbol, n::Int)

Takes a symbol and an integer and returns a vector of symbols with increasing numbers appended to the end. E.g., (:x, 3) -> [:x1, :x2, :x3]

The short-hand syntax s^n is also available, e.g., :x^3 == expand_symbol(:x, 3).

Useful to create signal names for named systems.

RobustAndOptimalControl.extended_controllerFunction
extended_controller(l::LQGProblem, L = lqr(l), K = kalman(l))

Returns an expression for the controller that is obtained when state-feedback u = -L(xᵣ-x̂) is combined with a Kalman filter with gain K that produces state estimates x̂. The controller is an instance of ExtendedStateSpace where C2 = -L, D21 = L and B2 = K.

The returned system has inputs [xᵣ; y] and outputs the control signal u. If a reference model R is used to generate state references xᵣ, the controller from e = ry - y -> u is given by

Ce = extended_controller(l)
Ce = named_ss(Ce; x = :xC, y = :u, u = [R.y; :y^l.ny]) # Name the inputs of Ce the same as the outputs of `R`.
connect([R, Ce]; u1 = R.y, y1 = R.y, w1 = [:ry^l.ny, :y^l.ny], z1=[:u])

Since the negative part of the feedback is built into the returned system, we have C = observer_controller(l) Ce = extended_controller(l) system_mapping(Ce) == -C`

RobustAndOptimalControl.extended_controllerMethod
extended_controller(K::AbstractStateSpace)

Takes a controller and returns an ExtendedStateSpace version which has augmented input [r; y] and output y (z output is 0-dim).

RobustAndOptimalControl.find_lftMethod
l = find_lft(sys::StateSpace{<:Any, <:StaticParticles{<:Any, N}}, δ) where N

Given an systems sys with uncertain coefficients in the form of StaticParticles, find a lower linear fractional transformation M such that lft(M, δ) ≈ sys.

δ can be either the source of uncertainty in sys, i.e., a vector of the unique uncertain parameters that were used to create sys. These should be constructed as uniform randomly distributed particles for most robust-control theory to be applicable. δ can also be an integer, in which case a numer of δ sources of uncertainty are automatically created. This could be used for order reduction if the number of uncertainty sources in sys is large.

Note, uncertainty in sys is only supported in A and B, C and D must be deterministic.

Returns l::LFT that internaly contains all four blocks of M as well as δ. Call ss(l,sys) do obtain lft(M, δ) ≈ sys.

Call Matrix(l) to obtain M = [M11 M12; M21 M22]

RobustAndOptimalControl.glover_mcfarlaneFunction
K, γmin = glover_mcfarlane(G::AbstractStateSpace{Continuous}, γ = 1.1)

Design a controller for G that maximizes the stability margin ϵ = 1/γ with normalized coprime factor uncertainty using the method of Glover and McFarlane

γ = 1/ϵ = ||[K;I] inv(I-GK)inv(M)||∞ G = inv(M + ΔM)*(N + ΔN)

γ is given as a relative factor above γmin and must be greater than 1, i.e., if γ = 1.1, the controller will be designed for γ = 1.1*γmin.

We want γmin ≥ 1 as small as possible, and we usually require that min is less than 4, corresponding to 25% allowed coprime uncertainty.

Performance modeling is incorporated in the design by calling glover_mcfarlane on the shaped system W2*G*W1 and then forming the controller as W1*K*W2. Using this formulation, traditional loop shaping can be done on W2*G*W1. Too many options? Select W2 diagonal with a small weight for less important outputs. Skogestad gives the following general advice:

  1. Scale the plant outputs and inputs. This is very important for most design procedures and is sometimes forgotten. In general, scaling improves the conditioning of the design problem, it enables meaningful analysis to be made of the robustness properties of the feedback system in the frequency domain, and for loop-shaping it can simplify the selection of weights. There are a variety of methods available including normalization with respect to the magnitude of the maximum or average value of the signal in question. Scaling with respect to maximum values is important if the controllability analysis of earlier chapters is to be used. However, if one is to go straight to a design the following variation has proved useful in practice: (a) The outputs are scaled such that equal magnitudes of cross-coupling into each of the outputs is equally undesirable. (b) Each input is scaled by a given percentage (say 10%) of its expected range of operation. That is, the inputs are scaled to reflect the relative actuator capabilities. An example of this type of scaling is given in the aero-engine case study of Chapter 12.

  2. Order the inputs and outputs so that the plant is as diagonal as possible. The relative gain array can be useful here. The purpose of this pseudo-diagonalization is to ease the design of the pre- and post-compensators which, for simplicity, will be chosen to be diagonal.

    Next, we discuss the selection of weights to obtain the shaped plant Gs = W2GW1 where W1 = Wp Wa Wg

  3. Select the elements of diagonal pre- and post-compensators Wp and W2 so that the singular values of W2 G Wp are desirable. This would normally mean high gain at low frequencies, roll-off rates of approximately 20 dB/decade (a slope of about 1) at the desired bandwidth(s), with higher rates at high frequencies. Some trial and error is involved here. W2 is usually chosen as a constant, reflecting the relative importance of the outputs to be controlled and the other measurements being fed back to the controller. For example, if there are feedback measurements of two outputs to be controlled and a velocity signal, then W 2 might be chosen to be diag[1, 1, 0.1], where 0.1 is in the velocity signal channel. W p contains the dynamic shaping. Integral action, for low frequency performance; phase-advance for reducing the roll-off rates at crossover; and phase-lag to increase the roll-off rates at high frequencies should all be placed in Wp if desired. The weights should be chosen so that no unstable hidden modes are created in Gs.

  4. Optional: Introduce an additional gain matrix Wg cascaded with Wa to provide control over actuator usage. Wg is diagonal and is adjusted so that actuator rate limits are not exceeded for reference demands and typical disturbances on the scaled plant outputs. This requires some trial and error.

  5. Robustly stabilize the shaped plant Gs = W2GW1 , where W1 = Wp Wa Wg , using glover_mcfarlane. First, the maximum stability margin ϵmax = 1/γmin is calculated. If the margin is too small, ϵmax < 0.25, then go back to step 5 and modify the weights. Otherwise, a γ-suboptimal controller is synthesized. There is usually no advantage to be gained by using the optimal controller. When ϵmax > 0.25 (respectively γmin < 4) the design is usually successful. In this case, at least 25% coprime factor uncertainty is allowed, and we also find that the shape of the open-loop singular values will not have changed much after robust stabilization. A small value of ϵmax indicates that the chosen singular value loop-shapes are incompatible with robust stability requirements. That the loop-shapes do not change much following robust stabilization if γ is small (ϵ large), is justified theoretically in McFarlane and Glover (1990).

  6. Analyze the design and if all the specifications are not met make further modifications to the weights.

  7. Implement the controller. The configuration shown in below has been found useful when compared with the conventional set up. This is because the references do not directly excite the dynamics of Ks, which can result in large amounts of overshoot (classical derivative kick). The constant prefilter ensures a steady-state gain of 1 between r and y, assuming integral action in W1 or G (note, the K returned by this function has opposite sign compared to that of Skogestad, so we use negative feedback here).

       ┌─────────┐      ┌────────┐      ┌────────┐
    r  │         │ -  us│        │  u   │        │  y
   ───►│(K*W2)(0)├──+──►│   W1   ├─────►│   G    ├────┬──►
       │         │  │+  │        │      │        │    │
       └─────────┘  │   └────────┘      └────────┘    │
                    │                                 │
                    │                                 │
                    │   ┌────────┐      ┌────────┐    │
                    │   │        │  ys  │        │    │
                    └───┤   K    │◄─────┤   W2   │◄───┘
                        │        │      │        │
                        └────────┘      └────────┘

Example:

Example 9.3 from the reference below.

using RobustAndOptimalControl, ControlSystems, Plots
G = tf(200, [10, 1])*tf(1, [0.05, 1])^2     |> ss
Gd = tf(100, [10, 1])                       |> ss
W1 = tf([1, 2], [1, 1e-6])                  |> ss
Gs = G*W1
Ks, γmin = glover_mcfarlane(Gs, 1.1)
@test γmin ≈ 2.34 atol=0.005

bodeplot([G, Gs, Gs*Ks]) |> display

plot( step(Gd*feedback(1, G*W1), 3))
plot!(step(Gd*feedback(1, G*W1*Ks), 3)) |> display

nyquistplot([G*W1, G*W1*Ks], ylims=(-2,1), xlims=(-2, 1), Ms_circles=1.5) |> display

Anti-windup can be added to W1 but putting W1 on Hanus form after the synthesis, see hanus.

Ref: Sec 9.4.1 of Skogestad, "Multivariable Feedback Control: Analysis and Design"

RobustAndOptimalControl.h2synthesizeFunction
K, Cl = h2synthesize(P::ExtendedStateSpace, γ = nothing)

Synthesize H₂-optimal controller K and calculate the closed-loop transfer function from w to z. Ref: Cha. 14.5 in Robust and Optimal Control.

If γ = nothing, use the formulas for H₂ in Ch 14.5. If γ is a large value, the H∞ formulas are used. As γ → ∞, these two are equivalent. The h∞ formulas do a coordinate transfromation that handles slightly more general systems so if you run into an error, it might be worth trying setting γ to something large, e.g., 1000.

RobustAndOptimalControl.hanusMethod
hanus(W)

Return Wh on Hanus form. Wh has twice the number of inputs, where the second half of the inputs are "actual inputs", e.g., potentially saturated. This is used to endow W with anti-windup protection. W must have an invertable D matrix and be minimum phase.

Ref: Sec 9.4.5 of Skogestad, "Multivariable Feedback Control: Analysis and Design"

RobustAndOptimalControl.hinfassumptionsMethod
flag = function hinfassumptions(P::ExtendedStateSpace; verbose=true)

Check the assumptions for using the γ-iteration synthesis in Theorem 1. In future revisions, we could suggest possible changes to P should the system not be feasible for synthesis. However, this has not been too successful so far..!

RobustAndOptimalControl.hinfpartitionMethod
P = hinfpartition(G, WS, WU, WT)

Transform a SISO or MIMO system G, with weighting functions WS, WU, WT into and LFT with an isolated controller, and write the resulting system, P(s), on a state-space form. Valid inputs for G are transfer functions (with dynamics, can be both MIMO and SISO, both in tf and ss forms). Valid inputs for the weighting functions are empty arrays, numbers (static gains), and LTISystems.

RobustAndOptimalControl.hinfsignalsMethod
hinfsignals(P::ExtendedStateSpace, G::LTISystem, C::LTISystem)

Use the extended state-space model, a plant and the found controller to extract the closed loop transfer functions operating solely on the state-space.

  • Pcl : w → z : From input to the weighted functions
  • S : w → e : From input to error
  • CS : w → u : From input to control
  • T : w → y : From input to output
RobustAndOptimalControl.hinfsynthesizeMethod
flag, K, γ, mats = hinfsynthesize(P::ExtendedStateSpace; maxIter=20, interval=(2/3,20), verbose=true)

Computes an H-infinity optimal controller K for an extended plant P such that ||F_l(P, K)||∞ < γ for the largest possible γ given P. The routine is known as the γ-iteration, and is based on the paper "State-space formulae for all stabilizing controllers that satisfy an H∞-norm bound and relations to risk sensitivity" by Glover and Doyle. See the Bib-entry below [1] above.

RobustAndOptimalControl.hsvdMethod
hsvd(sys::AbstractStateSpace{Continuous})

Return the Hankel singular values of sys, computed as the eigenvalues of QP Where Q and P are the Gramians of sys.

RobustAndOptimalControl.input_comp_sensitivityMethod
input_comp_sensitivity(P,C)
input_comp_sensitivity(l::LQGProblem)

Transfer function from load disturbance to control signal.

  • "Input" signifies that the transfer function is from the input of the plant.
  • "Complimentary" signifies that the transfer function is to an output (in this case controller output)
         ▲
         │e₁
         │  ┌─────┐
d₁────+──┴──►  P  ├─────┬──►e₄
      │-    └─────┘     │
      │                 │
      │     ┌─────┐     │
 e₂◄──┴─────┤  C  ◄──┬──+───d₂
            └─────┘  │
                     │e₃
                     ▼
  • Input sensitivity is the transfer function from d₁ to e₁, (I + CP)⁻¹
  • Output sensitivity is the transfer function from d₂ to e₂, (I + PC)⁻¹
  • Input complementary sensitivity is the transfer function from d₁ to e₂, CP(I + CP)⁻¹
  • Output complementary sensitivity is the transfer function from d₂ to e₄, PC(I + PC)⁻¹
RobustAndOptimalControl.input_sensitivityMethod
input_sensitivity(P, C)
input_sensitivity(l::LQGProblem)

Transfer function from load disturbance to total plant input.

  • "Input" signifies that the transfer function is from the input of the plant.
         ▲
         │e₁
         │  ┌─────┐
d₁────+──┴──►  P  ├─────┬──►e₄
      │-    └─────┘     │
      │                 │
      │     ┌─────┐     │
 e₂◄──┴─────┤  C  ◄──┬──+───d₂
            └─────┘  │
                     │e₃
                     ▼
  • Input sensitivity is the transfer function from d₁ to e₁, (I + CP)⁻¹
  • Output sensitivity is the transfer function from d₂ to e₂, (I + PC)⁻¹
  • Input complementary sensitivity is the transfer function from d₁ to e₂, CP(I + CP)⁻¹
  • Output complementary sensitivity is the transfer function from d₂ to e₄, PC(I + PC)⁻¹
RobustAndOptimalControl.named_ssMethod
named_ss(sys::AbstractStateSpace{T}; x, u, y)

Create a NamedStateSpace system. This kind of system uses names rather than integer indices to refer to states, inputs and outputs

Arguments:

  • sys: A system to add names to.
  • x: A list of symbols with names of the states.
  • u: A list of symbols with names of the inputs.
  • y: A list of symbols with names of the outputs.

Default names of signals if none are provided are x,u,y.

Example

```julia G1 = ss(1,1,1,0) G2 = ss(1,1,1,0) s1 = namedss(G1, x = :x, u = :u1, y=:y1) s2 = namedss(G2, x = :z, u = :u2, y=:y2)

s1[:y1, :u1] # Index using symbols

fb = feedback(s1, s2, r = :r) # ````

RobustAndOptimalControl.output_comp_sensitivityMethod
output_comp_sensitivity(P,C)
output_comp_sensitivity(l::LQGProblem)

Transfer function from measurement noise / reference to plant output.

  • "output" signifies that the transfer function is from the output of the plant.
  • "Complimentary" signifies that the transfer function is to an output (in this case plant output)
         ▲
         │e₁
         │  ┌─────┐
d₁────+──┴──►  P  ├─────┬──►e₄
      │-    └─────┘     │
      │                 │
      │     ┌─────┐     │
 e₂◄──┴─────┤  C  ◄──┬──+───d₂
            └─────┘  │
                     │e₃
                     ▼
  • Input sensitivity is the transfer function from d₁ to e₁, (I + CP)⁻¹
  • Output sensitivity is the transfer function from d₂ to e₂, (I + PC)⁻¹
  • Input complementary sensitivity is the transfer function from d₁ to e₂, CP(I + CP)⁻¹
  • Output complementary sensitivity is the transfer function from d₂ to e₄, PC(I + PC)⁻¹
RobustAndOptimalControl.output_sensitivityMethod
output_sensitivity(P, C)
output_sensitivity(l::LQGProblem)

Transfer function from measurement noise / reference to control signal.

  • "output" signifies that the transfer function is from the output of the plant.
         ▲
         │e₁
         │  ┌─────┐
d₁────+──┴──►  P  ├─────┬──►e₄
      │-    └─────┘     │
      │                 │
      │     ┌─────┐     │
 e₂◄──┴─────┤  C  ◄──┬──+───d₂
            └─────┘  │
                     │e₃
                     ▼
  • Input sensitivity is the transfer function from d₁ to e₁, (I + CP)⁻¹
  • Output sensitivity is the transfer function from d₂ to e₂, (I + PC)⁻¹
  • Input complementary sensitivity is the transfer function from d₁ to e₂, CP(I + CP)⁻¹
  • Output complementary sensitivity is the transfer function from d₂ to e₄, PC(I + PC)⁻¹
RobustAndOptimalControl.sensitivityMethod

See output_sensitivity

         ▲
         │e₁
         │  ┌─────┐
d₁────+──┴──►  P  ├─────┬──►e₄
      │-    └─────┘     │
      │                 │
      │     ┌─────┐     │
 e₂◄──┴─────┤  C  ◄──┬──+───d₂
            └─────┘  │
                     │e₃
                     ▼
  • Input sensitivity is the transfer function from d₁ to e₁, (I + CP)⁻¹
  • Output sensitivity is the transfer function from d₂ to e₂, (I + PC)⁻¹
  • Input complementary sensitivity is the transfer function from d₁ to e₂, CP(I + CP)⁻¹
  • Output complementary sensitivity is the transfer function from d₂ to e₄, PC(I + PC)⁻¹
RobustAndOptimalControl.specificationplotFunction
specificationplot([S,CS,T], [WS,WU,WT])

This function visualizes the control synthesis using the hInf_synthesize with the three weighting functions {WS(jω), WU(jω), WT(jω)} inverted and scaled by γ, against the corresponding transfer fucntions {S(jω), C(jω)S(jω), T(jω)}, to verify visually that the specifications are met. This may be run using both MIMO and SISO systems.

RobustAndOptimalControl.sumblockMethod
sumblock(ex::String; Ts = 0, n = 1)

Create a summation node that sums (or subtracts) vectors of length n.

Arguments:

  • Ts: Sample time
  • n: The length of the input and output vectors. Set n=1 for scalars.

Examples:

julia> sumblock("uP = vf + yL")
NamedStateSpace{Continuous, Int64}
D = 
 1  1

With state  names: 
     input  names: vf yL
     output names: uP


julia> sumblock("x_diff = xr - xh"; n=3)
NamedStateSpace{Continuous, Int64}
D = 
 1  0  0  -1   0   0
 0  1  0   0  -1   0
 0  0  1   0   0  -1

With state  names: 
     input  names: xr1 xr2 xr3 xh1 xh2 xh3
     output names: x_diff1 x_diff2 x_diff3
     

julia> sumblock("a = b + c - d")
NamedStateSpace{Continuous, Int64}
D = 
 1  1  -1

With state  names: 
     input  names: b c d
     output names: a