Feedback Linearization

Feedback linearization is a nonlinear control technique that transforms a nonlinear system into an equivalent linear system through a change of variables and a suitable control input. This technique is particularly useful for controlling nonlinear systems where the nonlinearities can be exactly canceled by the control input.

Theory

Consider a nonlinear system $P$ in input-affine form:

\[\begin{aligned} \dot{x} &= f(x) + g(x)u \\ y &= h(x) \end{aligned}\]

where $x \in \mathbb{R}^{n_x}$ is the state vector, $u \in \mathbb{R}^{n_u}$ is the control input, $y \in \mathbb{R}^{n_y = n_u}$ is the output, and $f(x)$, $g(x)$, and $h(x)$ are smooth vector fields.

The goal of feedback linearization is to find a control law that renders the input-output behavior of the system linear. This is achieved through:

  1. Computing Lie derivatives to determine the relative degree $\nu$
  2. Constructing a nonlinear coordinate transformation$z = T(x)$ where $z = [y, \dot{y}, \ddot{y}, ..., y^{(\nu-1)}]$
  3. Applying a nonlinear control law that cancels the nonlinearities

Relative Degree

The relative degree $\nu$ is the number of times the output $y$ must be differentiated before the input $u$ appears explicitly. It is determined by:

  • $L_g h(x) = 0$ (input doesn't appear after first differentiation)
  • $L_g L_f h(x) = 0$ (input doesn't appear after second differentiation)
  • ...
  • $L_g L_f^{\nu-1} h(x) \neq 0$ (input appears after $\nu$ differentiations)

where $L_f h$ denotes the Lie derivative of $h$ along the vector field $f$.

Feedback Linearizing Control Law

The feedback linearizing control law is given by:

\[u = (L_g L_f^{\nu-1} h)^{-1}(v - L_f^\nu h)\]

where $v$ is the new input for the linearized system, and the system from $v$ to $z$ behaves as a chain of $\nu$ integrators. A standard linear controller, such as a state-feedback controller designed using pole placement, may then be used to control the apparent linear system. Finally, the output of the linear controller is transformed using the nonlinear transformation above to produce the true control input $u$ to the nonlinear plant $P$.

Control Architecture

The complete control architecture can be visualized as:

r    ┌─────┐     ┌─────┐      ┌─────┐      
────►│     │ v   │     │  u   │     │ x    
     │  Cₗ ├────►│ Cₙₗ ├─────►│  P  ├─┬───►
  ┌─►│     │     │     │      │     │ │    
  │  └─────┘     └─────┘      └─────┘ │    
  │                 ▲                 │    
  │        ┌─────┐  │                 │    
  │      z │     │  │                 │    
  └────────┤  T  │◄─┴─────────────────┘    
           │     │                         
           └─────┘                         

Where:

  • $C_l$: Linear controller (e.g., PID, state feedback)
  • $C_{nl}$: Nonlinear transformation (feedback linearization)
  • $P$: Original nonlinear plant
  • $T$: Nonlinear coordinate transformation $z = T(x)$

Implementation in DyadControlSystems

DyadControlSystems provides the following functions for feedback linearization:

Examples

Example 1: Simple Third-Order System with PDD² Controller

This example demonstrates feedback linearization of a third-order system with sinusoidal nonlinearity, implementing a PDD² (proportional-derivative-derivative²) controller for trajectory tracking. The reference trajectory is a sinusoidal function of time.

using DyadControlSystems
using Symbolics
using Plots

@variables x[1:3] u
x = collect(x)
f = [
    x[2]
    sin(x[3])
    0
]
g = [0; 0; 1]
h = x[1]

result = feedbacklin(f, g, h, x)
FeedbackLinearizationResult:
  Relative degree: 3
  
  Lie derivatives:
    L_f^0 h = h = x[1]
    L_f^1 h = x[2]
    L_f^2 h = sin(x[3])
  
  Control-related terms:
    L_f^ν h = 0
    L_g L_f^(ν-1) h = cos(x[3])
  
  Feedback linearizing control law:
    u = v / cos(x[3])
  
  New input variable: v

The result is a FeedbackLinearizationResult object containing all the necessary information for simulating the feedback linearized system, but in symbolic form. To generate code for simulation, we call build_function. We also implement the linear state-feedback controller in a function v_func. This computes the linear control output v given the nonlinear state x, we call the nonlinear coordinate transformation $T$ inside this function to perform the change of coordinates from $x$ to $z$. Since $z$ is comprised of the output $y$ and it's derivatives, we may implement a PD-style control law as state feedback in this linear coordinate system. The function simulate automatically transforms $v$ into $u$, otherwise one would call u = sys.u(x, v, p, t) to do this.

sys = Symbolics.build_function(result)

# Initial condition
x0 = [0.5, 0.0, π/4]

# Reference input: sinusoidal trajectory tracking with PDD² feedback
v_func = function (x,p,t)
    y, yd, ydd = sys.T(x, p, t) # Perform coordinate change from x to z = [y, ẏ, ÿ]
    e = sin(0.5*t) - y
    ed = 0.5*cos(0.5*t) - yd
    edd = -0.5^2 * sin(0.5*t) - ydd
    5*(e + 1*ed + 1*edd) # proportional, derivative, and second derivative feedback
end

# Simulate
sol = DyadControlSystems.simulate(sys, x0, v_func, (0.0, 30.0); Ts=0.01)
plot(sol, ploty=false, plotx=false, plotv=true, plotu=true, plotz=true,
     size=(800, 600), layout=(5,1))
plot!(sol.t, sin.(0.5.*sol.t), lab="reference", sp=3, l=(:black, :dash))
plot!(sol.t, @.(0.5cos(0.5sol.t)), lab="reference", sp=4, l=(:black, :dash))
plot!(sol.t, @.(-0.5^2 * sin(0.5*sol.t)), lab="reference", sp=5, l=(:black, :dash))
Example block output

The system has relative degree 3, meaning the output must be differentiated three times before the input appears. The PDD² controller provides feedback on the tracking error and its first two derivatives, effectively controlling all three linearized state variables.

In this example, we provided reference trajectories for the linear coordinates $z$. The linear coordinates correspond to the output $y = h(x)$ and its derivatives.

Example 2: Pole Placement

This example shows how to use pole placement for the linearized system to achieve desired closed-loop dynamics.

using LinearAlgebra

@variables x[1:3] u
x = collect(x)
f = [
    sin(x[2]) + (x[2]+1)*x[3]
    x[1]^5 + x[3]
    x[1]^2
]
g = [0; 0; 1]
h = x[1]

result = feedbacklin(f, g, h, x)
sys = Symbolics.build_function(result; simplify=false)
L = place(sys, [-1 + 0.1im, -1 - 0.1im]) # Place poles for the linearized system, L is feedback gain matrix

# Initial condition
x0 = [0.1, 0.1, 0.1]

# State feedback control using pole placement
v_func = function (x, p, t)
    -(L*sys.T(x, p, t))[]
end

# Simulate
sol = DyadControlSystems.simulate(sys, x0, v_func, (0.0, 15.0))
plot(sol, ploty=false, plotx=true, plotv=true, plotu=false, plotz=false,
     size=(800, 800), framestyle=:zerolines)
Example block output

The pole placement technique designs a state feedback controller for the linearized system. The poles are placed at $-1 \pm 0.1i$, providing a well-damped response.

Example 3: MIMO System

Feedback linearization can also be applied to Multi-Input Multi-Output (MIMO) systems, provided the number of inputs is at least equal to the number of outputs.

@variables x[1:4] u[1:2]
x = collect(x)
u = collect(u)
p = []

# System with two inputs affecting different parts
f = [
    x[2]
    -x[1] + x[3]^2
    x[4]
    -x[3] + sin(x[1])
]

# Input matrix: each input affects different states
g = [
    0 0
    1 0
    0 0
    0 1
]

# Two outputs
h = [x[1], x[3]]
result = feedbacklin(f, g, h, x)
sys = Symbolics.build_function(result)
x0 = ones(4) # Initial condition

# Reference tracking for both outputs
v_func = function (x,p,t)
    y = sys.T(x, p, t)
    v1 = [2,0] # reference for output 1 and its derivative
    v2 = [3,0] # reference for output 2 and its derivative
    V = [v1 v2]
    return sum(V .- y, dims=1) |> vec
end

# Simulate
sol = DyadControlSystems.simulate(sys, x0, v_func, (0.0, 20.0))
plot(sol, ploty=true, plotx=false, plotv=true, plotu=true, plotz=false,
     size=(800, 800))
hline!([2 3], lab="reference", sp=(1:2)', c=:black, ls=:dash)
Example block output

In this MIMO example, the system has two inputs and two outputs. Each output is feedback linearized independently, and the resulting controller tracks different reference values for each output.

Example 4: ModelingToolkit Integration

DyadControlSystems seamlessly integrates with ModelingToolkit.jl for physical modeling. This example demonstrates feedback linearization on a system defined using ModelingToolkit, recreating the third-order system from Example 1.

using ModelingToolkit
using Symbolics
using ModelingToolkit: t_nounits as t, D_nounits as D

# Define the same third-order system using ModelingToolkit
@mtkmodel ThirdOrderSystem begin
    @variables begin
        x1(t) = 0.5
        x2(t) = 0.0
        x3(t) = π/4
        u(t)         # control input
        y(t)         # output
    end
    @equations begin
        D(x1) ~ x2
        D(x2) ~ sin(x3)
        D(x3) ~ u
        y ~ x1
    end
end

# Create and compile the model
@named model = ThirdOrderSystem()
cmodel = complete(model)
inputs = [cmodel.u]
h_mtk = cmodel.y # Specify controlled output, this currently must be a single variable

result_mtk, x, p, sys = feedbacklin(model, inputs, h_mtk)

sys_mtk = build_function(result_mtk; p)

# Use the same PDD² controller from Example 1
v_func_mtk = function (x,p,t)
    y, yd, ydd = sys_mtk.T(x, p, t)
    e = sin(0.5*t) - y
    ed = 0.5*cos(0.5*t) - yd
    edd = -0.5^2 * sin(0.5*t) - ydd
    5*(e + 1*ed + 1*edd)
end

# Simulate with same initial condition
x0_mtk = ModelingToolkit.get_u0(sys, []) # MTK uses u to denote state in this context
sol_mtk = DyadControlSystems.simulate(sys_mtk, x0_mtk, v_func_mtk, (0.0, 30.0))

# Plot and compare with symbolic version
plot(sol_mtk, ploty=true, plotx=false, plotv=false, plotu=true, plotz=false,
     size=(800, 400), layout=(2,1))
plot!(sol_mtk.t, sin.(0.5.*sol_mtk.t), lab="reference", sp=1, ls=:dash, c=:black)
Example block output

ModelingToolkit offers several advantages for control system design:

  • Physical modeling: Define systems using physical equations directly
  • Automatic differentiation: Symbolic manipulation and automatic Jacobian computation
  • Parameter handling: Easy parameter management and sensitivity analysis
  • Component-based modeling: Build complex systems from reusable components

The input_affine_form function extracts the control-affine representation $\dot{x} = f(x) + g(x)u$ from the ModelingToolkit equations, making it compatible with the feedback linearization procedure.

Practical Considerations

When to Use Feedback Linearization

Feedback linearization is most effective when:

  • The system dynamics are known accurately
  • Accurate state estimation is available
  • The nonlinearities are smooth and differentiable

Limitations

  • Can be sensitive to model uncertainties
  • May cancel beneficial nonlinearities
  • Not applicable to all nonlinear systems (requires specific structure)

Index

DyadControlSystems.FeedbackLinearizationSolutionType
FeedbackLinearizationSolution

Solution structure containing the results of a feedback linearization simulation.

Fields:

  • sys: The FeedbackLinearizedSystem used for simulation
  • t: Time vector
  • x: State trajectories (nx × nt matrix)
  • y: Output trajectories (ny × nt matrix)
  • v: Control input v trajectories (nu × nt matrix)
  • u: Actual control u trajectories (computed from v, nu × nt matrix)
  • z: Linearized coordinates z = T(x) (nz × nt matrix)

Plot

The solution object can be plotted:

plot(sol; plotx=true, ploty=true, plotv=true, plotu=true, plotz=true)

where:

  • plotx: Plot the state trajectories
  • ploty: Plot the output trajectories
  • plotv: Plot the new control input v
  • plotu: Plot the actual control input u
  • plotz: Plot the linearized coordinates z
source
DyadControlSystems.FeedbackLinearizedSystemType
FeedbackLinearizedSystem

A structure containing compiled functions for simulating a feedback linearized system.

Fields:

  • f: Drift dynamics function f(x,p,t)
  • g: Input dynamics function g(x,p,t)
  • u: Feedback control law u(x,v,p,t)
  • T: Lie derivatives T(x,p,t), this function performs the nonlinear change of coordinates from x to z (y and its derivatives), i.e., the system from v to z behaves like a chain of integrators.
  • res: Original FeedbackLinearizationResult

Example:

The following implements a PDD² controller for a simple third-order system

@variables x[1:3] u
x = collect(x)
f = [
    x[2]
    sin(x[3])
    0
]
g = [0; 0; 1]
h = x[1]

result = feedbacklin(f, g, h, x)
sys = Symbolics.build_function(result; simplify=false)

# Initial condition
x0 = [0.5, 0.0, π/4]

# Reference input: sinusoidal
v_func = function (x,p,t)
    y, yd, ydd = sys.T(x, p, t) # Perform coordinate change from x to z = [y, ẏ, ÿ]
    e = sin(0.5*t) - y
    ed = 0.5*cos(0.5*t) - yd
    edd = -0.5^2 * sin(0.5*t) - ydd
    5*(e + 1*ed + 1edd) # proportional, derivative, and second derivative feedback
end

# Simulate
sol = DyadControlSystems.simulate(sys, x0, v_func, (0.0, 30.0); Ts=0.01)
plot(sol.t, sol.y')
plot!(sol.t, sin.(0.5.*sol.t), lab="ref")

Pole placement

To perform pole placement in the new coordinate system, call

L = place(sys, desired_poles)

This returns feedback gain L that can be multiplied by v = -L*sys.T(x, p, t) to compute the control signal for the linearized system. This, in turn, should be passed through sys.u(x, v, t) to compute the actual control signal u for the original nonlinear system.

Example

@variables x[1:3] u
x = collect(x)
f = [
    sin(x[2]) + (x[2]+1)*x[3]
    x[1]^5 + x[3]
    x[1]^2
]
g = [0; 0; 1]
h = x[1]

result = feedbacklin(f, g, h, x)
sys = Symbolics.build_function(result; simplify=false)
L = place(sys, [-0.5 + 0.1im, -0.5 - 0.1im]) # Place poles for the linearized system
# Initial condition
x0 = [0.1, 0.1, 0.1]

# Reference input: ramp function
v_func = function (x, p, t)
    -(L*sys.T(x, p, t))[]
end

# Simulate
sol = DyadControlSystems.simulate(sys, x0, v_func, (0.0, 15.0); Ts=0.01)
plot(sol.x')

MIMO system

For MIMO systems, the number of outputs and inputs must be identical, i.e., the length of h must be the width of g. The output of T is in this case a matrix where the width is the number of outputs

@variables x[1:4] u[1:2]
x = collect(x)
u = collect(u)
p = []

# System with two inputs affecting different parts
f = [
    x[2]
    -x[1] + x[3]^2
    x[4]
    -x[3] + sin(x[1])
]

# Input matrix: each input affects different states
g = [
    0 0
    1 0
    0 0
    0 1
]

# Test first output h1 = x[1] with full g matrix
h = [x[1], x[3]]
result = feedbacklin(f, g, h, x)
sys = build_function(result)
x0 = ones(4) # Initial condition

v_func = function (x,p,t)
    y = sys.T(x, p, t)
    v1 = [2,0] # reference for output 1 and its derivative
    v2 = [3,0]
    V = [v1 v2]
    return sum(V .- y, dims=1) |> vec
end

# Simulate
sol = DyadControlSystems.simulate(sys, x0, v_func, (0.0, 50.0); Ts=0.01)
plot(sol.y', layout=2)
hline!([2 3], lab="ref")
source
ControlSystemsBase.placeMethod
L = place(sys::FeedbackLinearizedSystem, desired_poles)

Pole placement for feedback-linearized systems.

For systems with multiple outputs, desired_poles must be a vector of vectors, one for each output. The inner vectors indicate the desired pole locations for each output system, the length of this is given by the relative degree. In this case, use the returned matrix like this

z = T(x) # Perform nonlinear coordinate change
v = -L*vec(z) # Compute linear control input
source
DyadControlSystems.feedbacklinMethod
result = feedbacklin(f, g, h, x)

Feedback linearize a SISO input-output dynamical system on input-affine form

\[ẋ = f(x) + g(x)u \ y = h(x)\]

The feedback-linearizing control law is given by

\[u = \dfrac{1}{L_g L_f^{ν-1} h} (v - L_f^ν h)\]

which is realized by

u = 1/LgLfh * (v - Lfh)

where v is the new input and Lfh and LgLfh are symbolic expressions in x representing the Lie-derivatives. The linearized system from v to z will behave like a chain of ν integrators.

Returns a FeedbackLinearizationResult containing:

  • Lfh: Lie derivative L_f^ν h
  • LgLfh: Lie derivative Lg Lf^{ν-1} h
  • T: All Lie derivatives from order 0 to ν-1. This is an array of expressions of x, the first element of which is $y$, followed by derivatives of $y$ in increasing order.
  • relative_degree: The relative degree ν
  • control_law: The feedback linearizing control law

The complete controller can be represented as the composition of the nonlinear coordinate change operator $z = [y, ẏ, ÿ, ...] = T(x)$, a linear controller that acts on the linearized state $z$ (such as a state-feedback controller v = L*z), and a nonlinear transformation of v to form u = Cₙₗ(x, v). The complete composition may thus look like

\[u = Cₙₗ(x, L*T(x))\]

r    ┌─────┐     ┌─────┐      ┌─────┐      
────►│     │ v   │     │  u   │     │ x    
     │  Cₗ ├────►│ Cₙₗ ├─────►│  P  ├─┬───►
  ┌─►│     │     │     │      │     │ │    
  │  └─────┘     └─────┘      └─────┘ │    
  │                 ▲                 │    
  │        ┌─────┐  │                 │    
  │      z │     │  │                 │    
  └────────┤  T  │◄─┴─────────────────┘    
           │     │                         
           └─────┘                         

Arguments:

  • f: A vector of symbolic expressions representing f(x)
  • g: An array of symbolic expressions representing g(x) (a static numeric array is okay if g is not a function of x)
  • h: A symbolic expression representing h(x)
  • x: A vector of the symbolic x variables
source
DyadControlSystems.feedbacklinMethod
fl_result, x_sym, p_sym, simplified_sys = feedbacklin(sys::System, inputs, h)

Perform feedback linearization on a ModelingToolkit.System.

Arguments:

  • sys: A ModelingToolkit.System representing the dynamical system
  • inputs: A vector of input variables (e.g., [u] or [u1, u2])
  • h: A symbolic expression representing the output variable (e.g., y or [y1, y2])

Returns:

  • fl_result: A FeedbackLinearizationResult containing the results of the feedback linearization
  • x_sym: A vector of symbolic state variables
  • p_sym: A symbolic parameter object
  • simplified_sys: The compiled ModelingToolkit.System used for the linearization
source
DyadControlSystems.reldegMethod
reldeg(f, g, h, x)

Compute the relative degree of the system defined by f, g, and h.

The relative degree is the number of times the output h must be differentiated before the input appears explicitly.

Arguments:

  • f: System drift dynamics
  • g: Input dynamics
  • h: Output function
  • x: State variables
source
DyadControlSystems.simulateMethod
simulate(sys::FeedbackLinearizedSystem, x0, v_func, tspan; Ts=0.01, p=nothing)

Simulate a feedback linearized control-affine system.

The system dynamics are:

\[ẋ = f(x) + g(x)u\]

where the control law u is computed from the feedback linearization.

Arguments:

  • sys: FeedbackLinearizedSystem from build_function
  • x0: Initial state vector
  • v_func: Reference input function v(t)
  • tspan: Time span (t0, tf) for simulation
  • Ts: Sample time (default 0.01)
  • p: Parameters (optional)

Returns:

  • sol: Solution object with fields t (time) and u (states)
source