Skip to main content

Nonlinear Systems and Conversion to LPV

Overview

LPVcore is able to represent nonlinear systems in state-space representation and to convert these to LPV state-space representations. On this page, we give a brief overview of this functionality. Firstly, we give an overview of the nonlinear state-space systems object in LPVcore. After that, we give an overview of the conversion methods that are available.

Nonlinear State-Space Systems

General information

Suggested example

An LPVcore.nlss object allows the user to represent a nonlinear dynamical system in state-space representation of the form

ξx=f(x,u),y=h(x,u),\begin{aligned} \xi x &= f(x,u), \\[4pt] y &= h(x,u), \end{aligned}

where xx is the state of the model, uu its input, and yy its output, ξ=ddt\xi=\frac{d}{dt} in continuous-time and ξ=q1\xi=q^{-1} in discrete-time. Moreover, ff and hh are nonlinear functions, represented by MATLAB function handles. For example, the following code snippet constructs a continuous-time nonlinear state-space representation:

f = @(x,u) [-x(2, :); 3.0 * (1 - x(1, :).^2) .* x(2, :) - x(1) + u];
h = @(x,u) x(1, :);

nx = 2;
nu = 1;
ny = 1;

sys = LPVcore.nlss(f, h, nx, nu, ny);

Discrete-time systems can be constructed by also providing a sampling time:

Ts = 0.1;
sys = LPVcore.nlss(f, h, nx, nu, ny, Ts);
Symbolic Representation

By default, it is assumed that the ff and hh functions of LPVcore.nlss objects can be represented symbolically, using the Symbolic Toolbox of MATLAB. This is done such that the Jacobian matrix functions of ff and hh can be computed, which are required for the conversion methods to LPV representations and/or compute the differential form of nonlinear system. However, in cases where this is not required and/or ff and hh cannot be represented symbolically, the internal symbolic representation can be disabled through the 7th argument of LPVcore.nlss. The following code snippet, creates an LPVcore.nlss object for which the internal symbolic representation is disabled:

sys = LPVcore.nlss(f, h, nx, nu, ny, Ts, false);
Simulink block

Nonlinear systems represented by LPVcore.nlss objects can also be used in Simulink by using the corresponding Nonlinear System block, located in LPVcore/Nonlinear in the Library Browser.

Available methods

The following methods are currently available for LPVcore.nlss objects:

MethodDescription
c2dDiscretization of continuous-time nonlinear state-space systems.
isctCheck if the LPVcore.nlss object is continuous-time.
simSimulation of the nonlinear state-space system.

Next, we give a more in-depth overview for some of these methods.

System discretization using c2d

Similar to the c2d command for linear state-space objects in the Control Systems Toolbox of MATLAB, the c2d function for LPVcore.nlss objects allows for the discretization of continuous-time nonlinear state-space objects to discrete-time. For example, the following code snippet constructs a nonlinear state-space representation and converts it to discrete-time:

% CT System
f = @(x,u) [-x(2, :); 3.0 * (1 - x(1, :).^2) .* x(2, :) - x(1) + u];
h = @(x,u) x(1, :);

nx = 2;
nu = 1;
ny = 1;

sys = LPVcore.nlss(f, h, nx, nu, ny);

% Convert to DT
Ts = 0.1;
sysdt = c2d(sys, Ts);

By default, the forward Euler method is used to convert the system to discrete-time. Currently there are two conversion methods available: Forward Euler, using the eul option, and fourth order Runge-Kutta method, using rk4. For both methods, the input signal is assumed to be constant between samples. The conversion method can be set by specifying a third argument in c2d. For example, the following code snippet converts a system to discrete-time based on a fourth order Runge-Kutta method:

Ts = 0.1;
sysdt = c2d(sys, Ts, 'rk4');

System simulation using sim

Nonlinear state-space objects can also be simulated for a given input trajectory and initial condition through the sim command. For example, the following code snippet creates a continuous-time nonlinear system an simulates it for a specified initial condition and input:

% System
f = @(x,u) [-x(2, :); 3.0 * (1 - x(1, :).^2) .* x(2, :) - x(1) + u];
h = @(x,u) x(1, :);

nx = 2;
nu = 1;
ny = 1;

sys = LPVcore.nlss(f, h, nx, nu, ny);

% Simulation
x0 = [0.5; 0];
u = @(t) sin(2 * pi * 0.5 * t);
tspan = [0, 10];

[y, t, x] = sim(sys, u, tspan, x0);

The third argument, specifies either the time-span of the simulation, or a vector of time instances. In the latter case, the output time values will be equal to the values of the input time vector. For example:

% Simulation
x0 = [0.5; 0];
u = @(t) sin(2 * pi * 0.5 * t);
tvec = linspace(0, 10, 1000)';

[y, t, x] = sim(sys, u, tvec, x0); % t == tvec

For a time vector of length N, the input to the system can also be given in terms of a vector of size [N, Nu], where Nu is the input size. In this case, it is assumed that the input is constant between time instances. See also the following code snippet:

% Simulation
x0 = [0.5; 0];
tvec = linspace(0, 10, 1000)';
u = sin(2 * pi * 0.5 * tvec);

[y, t, x] = sim(sys, u, tvec, x0);

For continuous-time systems, ode45 is used by default in order to simulate the response of the system. For continuous-time systems, other solvers, such as ode89 or ode15s to name a few, can also be specified through the 4th argument of the sim command. For example, see the following code snippet to simulate a continuous-time nonlinear system using ode89:

[y, t, x] = sim(sys, u, tvec, x0, 'ode89');

The differential form

The differential form (see also Section 5.2 in [1]) of a nonlinear system represents the dynamics of the variations along its trajectories. For a nonlinear system in state-space form, the dynamics of the differential form are given by:

ξ(δx)=A(x,u)δx+B(x,u)δu,δy=C(x,u)δx+D(x,u)δu,\begin{aligned} \xi (\delta x) &= A(x,u) \delta x + B(x,u) \delta u, \\[4pt] \delta y &= C(x,u) \delta x + D(x,u) \delta u, \end{aligned}

where δx\delta x, δu\delta u, and δy\delta y are the differential state, input, and output, respectively. Furthermore, A=fxA = \frac{\partial f}{\partial x} , B=fuB = \frac{\partial f}{\partial u} , C=hxC = \frac{\partial h}{\partial x} , and D=huD = \frac{\partial h}{\partial u}. The differential form of an LPVcore.nlss representation can be computed through the LPVcore.difform command, for example:

f = @(x,u) [-x(2, :); 3.0 * (1 - x(1, :).^2) .* x(2, :) - x(1) + u];
h = @(x,u) x(1, :);

nx = 2;
nu = 1;
ny = 1;

sys = LPVcore.nlss(f, h, nx, nu, ny);
sysdif = LPVcore.difform(sys);

Currently, the differential form is only used (internally) for computation of the Jacobians, AA, BB, CC, DD, which are used for the nonlinear to LPV system conversion methods. Nevertheless, a few methods are available for LPVcore.difform objects in order to check the dependency of the differential form on (which part of) xx and uu, using the dependency method, and to simplify the differential form, using the simplify method. We refer the reader to the help of these methods (help LPVcore.difform/dependency and help LPVcore.difform/simplify) and the examples in nonlinear/example_difform.m for more details.

Nonlinear to LPV System Conversion

Overview

The table below shows the current list of functions which can be used to convert an LPVcore.nlss object to an LPV system representation:

FunctionDescription
LPVcore.nlss2lpvssConvert LPVcore.nlss object to an LPVcore.lpvss object.
LPVcore.nlss2lpvgridssConvert LPVcore.nlss object to an LPVcore.lpvgridss object.

The conversion for both of these functions from a nonlinear state-space system to an LPV representation is based on the results in Theorem C.6.1 in [1]. This approach is based on the fundamental theorem of calculus, which uses an integral of the Jacobians of the functions ff and hh to factorize them.

Zero Solution

For the conversion methods, it is assumed that f(0,0)=0f(0,0) = 0 and h(0,0)=0h(0, 0) = 0. This might not always be the case, however, one can always compute offsets for the system such that this is the case. For example, assume that f(0,0)=xf(0,0) = x^* and h(0,0)=yh(0,0) = y^*. Then define f~(x,u)=f(x,u)x\tilde f(x,u) = f(x,u) - x^* and h~(x,u)=h(x,u)y.\tilde h(x,u) = h(x,u) - y^*. The system defined by f~\tilde f and h~\tilde h will then satisfy that f~(0,0)=0\tilde f(0,0) = 0 and h~(0,0)=0\tilde h(0,0) = 0.

Conversion using LPVcore.nlss2lpvss

The LPVcore.nlss2lpvss allows us to convert an LPVcore.nlss object to an LPVcore.lpvss object. For example, see the following code snippet:

f = @(x,u) [-x(2, :); 3.0 * (1 - x(1, :).^2) .* x(2, :) - x(1) + u];
h = @(x,u) x(1, :);

nx = 2;
nu = 1;
ny = 1;

sys = LPVcore.nlss(f, h, nx, nu, ny);

[sysLPV, map] = LPVcore.nlss2lpvss(sys);

The function returns both the LPV system as LPVcore.lpvss object and a scheduling-map, as LPVcore.schedmap object, which is a function mapping the states and inputs of the system to the scheduling-variable of the LPV representation, i.e.:

p=ψ(x,u),p = \psi(x,u),

where ψ\psi is the scheduling-map.

As aforementioned, the conversion is based on the fundamental theorem of calculus, which uses an integral of the Jacobians of the functions ff and hh. By default, this integral is computed numerically and contained in the scheduling-map. However, through the second argument of the LPVcore.nlss2lpvss function, the user can specify to compute this integral analytically, see for example the following code snippet:

f = @(x,u) [-x(2, :); 3.0 * (1 - x(1, :).^2) .* x(2, :) - x(1) + u];
h = @(x,u) x(1, :);

nx = 2;
nu = 1;
ny = 1;

sys = LPVcore.nlss(f, h, nx, nu, ny);

[sysLPV, map] = [sysLPV, map] = LPVcore.nlss2lpvss(sys, "analytical");
Analytical Option

Note that the "analytical" option for LPVcore.nlss2lpvss might not always work. Namely, dependending on the functions ff and hh, analytically computing an integral over the Jacobians might not always be possible.

Finally, remember that the LPV system resulting from the conversion is an LPVcore.lpvss object. This means that its AA, BB, CC, DD matrices are given by pmatrix objects of the form1:

M(p)=i=1NMiϕi(p).M(p) = \sum_{i=1}^{N} M_i \phi_i(p).

By default, due to the conversion procedure, the MiM_i matrices are binary matrices and ϕi(p)=pi\phi_i(p) = p_i. This means that all the functional dependency is 'hidden' inside the scheduling-map ψ(x,u)\psi(x,u). Through the third argument of the LPVcore.nlss2lpvss function, the user can specify to use a more advanced option, which tries to factorize the dependency further, using "factor". For example,

f = @(x,u) [-x(2, :); 3.0 * (1 - x(1, :).^2) .* x(2, :) - x(1) + u];
h = @(x,u) x(1, :);

nx = 2;
nu = 1;
ny = 1;

sys = LPVcore.nlss(f, h, nx, nu, ny);
[sysLPV, map] = LPVcore.nlss2lpvss(sys, "analytical", "factor");

results in:

>> sysLPV.A

ans =

Continuous-time scheduling dependent 2x2 matrix

A1 + A2*(p1(t)) + A3*(p2(t))

A1 =
0 -1
-1 3

A2 =
0 0
-2 0

A3 =
0 0
0 -1

>> map.map

ans =

function_handle with value:

@(x,u)[x(1,:).*x(2,:);x(1,:).^2]

while

[sysLPV, map] = LPVcore.nlss2lpvss(sys, "analytical", "element");   % default for third argument

results in:

>> sysLPV.A

ans =

Continuous-time scheduling dependent 2x2 matrix

A1 + A2*(p1(t)) + A3*(p2(t))

A1 =
0 -1
0 0

A2 =
0 0
1 0

A3 =
0 0
0 1

>> map.map

ans =

function_handle with value:

@(x,u)[x(1,:).*x(2,:).*-2.0-1.0;-x(1,:).^2+3.0]

Notice that while both LPV representations are equivalent representation of the corresponding nonlinear state-space system, they have different AiA_i matrices and scheduling-maps ψ(x,u)\psi(x,u).

Conversion using LPVcore.nlss2lpvgridss

The LPVcore.nlss2lpvgridss function allows us to convert an LPVcore.nlss object to an LPVcore.lpvgridss object. For example, see the following code snippet:

f = @(x,u) [-x(2, :); 3.0 * (1 - x(1, :).^2) .* x(2, :) - x(1) + u];
h = @(x,u) x(1, :);

nx = 2;
nu = 1;
ny = 1;

sys = LPVcore.nlss(f, h, nx, nu, ny);

grid = struct('x1', [-1:0.5:1], 'x2', [-3:1:3], 'u', [-5:1:5]);
[sysLPVgrid] = LPVcore.nlss2lpvgridss(sys, grid);

The second argument of LPVcore.lpvgridss specifies the state and input grid-points for the gridded LPV representation. As shown in the example above, this grid is represented by a struct, which contains as a field-value pair for each state (two in the example) and for each input (one in the example). Each field for the state, should be denoted by xi, where i = 1, 2, ..., Nx corresponds to the i'th state, and similarly for the input. If there's only one state (or input), one can use x (or u in the input case).

Choosing your grid

When choosing your grid points, it is not always necessary to give multiple grid-points for each state and/or input to exactly represent the nonlinear system. For example, the nonlinear system in the examples above is linear in its input, that means that it's actually not necessary to grid over the input (as the LPV model will not vary in uu). Therefore, it's only necessary to grid over the states and inputs that appear nonlinearly in ff and hh. The states and inputs which appear nonlinearly can also be checked with the help of the differential form of the system and the dependency command (as the Jacobians of ff and hh will depend on states and inputs that appear nonlinearly). For example, for the system above, we have that:

>> dfsys = LPVcore.difform(sys);
>> dependency(dfsys)

State-space matrix dependency of differential form:

A: ( x1, x2 )
B: None
C: None
D: None

Based on this, we can see that we only have nonlinearities depending on x1x_1 and x2x_2. As we have to provide grid-points for all states and inputs for LPVcore.nlss2lpvgridss, we can achieve gridding over only x1x_1 and x2x_2 by specifying only one (arbitrary) grid-point for uu:

grid = struct('x1', [-1:0.5:1], 'x2', [-3:1:3], 'u', 0);
[sysLPVgrid] = LPVcore.nlss2lpvgridss(sys, grid);

Differential form conversion methods

It's also possible to convert differential forms of a nonlinear system to an LPV representation. An overview of the available functions for this is given in the following table:

FunctionDescription
LPVcore.difform2lpvssConvert LPVcore.difform object to an LPVcore.lpvss object.
LPVcore.difform2lpvgridssConvert LPVcore.difform object to an LPVcore.lpvgridss object.

These functions work similar to their counterpart LPVcore.nlss2{*} functions, therefore, we do not elaborate on these further.

Scheduling-map objects

As briefly mentioned, besides an LPV representation as an LPVcore.lpvss object, the LPVcore.nlss2lpvss and LPVcore.difform2lpvss functions also return a LPVcore.schedmap objects, which represents a scheduling-map of the form

p=ψ(x,u),p = \psi(x,u),

where ψ\psi is the scheduling-map.

A scheduling-map object can also be constructed manually by providing a function handle, corresponding to the map, and providing the corresponding number of scheduling-variables. For example, to represent the scheduling-map:

ψ(x,u)=[sin(x1)x1x2u2]\psi(x,u) = \begin{bmatrix} \sin(x_1)\\ x_1 x_2\\ u_2\end{bmatrix}

one would call:

schMap = LPVcore.schedmap(@(x,u) [sin(x(1, :)); x(1, :) .* x(2, :); u(2, :)], 3);

For this example, the map can be evaluated by manually calling schMap.map(x,u) or eval(schMap, x, u), where x and u are matrices of size [Nx, N] and [Nu, N], respectively. See for example, the following code snippet:

nx = 2;
nu = 2;
N = 10;

x = randn(nx, N);
u = randn(nu, N);

pTraj1 = schMap.map(x,u);
pTraj2 = eval(schMap, x, u); % pTraj1 == pTraj2

References

[1]: Koelewijn, P.J.W., 'Analysis and Control of Nonlinear Systems with Stability and Performance Guarantees: A Linear Parameter-Varying Approach', PhD Thesis, Electrical Engineering, Eindhoven, 2023.


  1. For the conversion methods, the extended scheduling signal is alway equal to the scheduling signal itself, i.e., ρ(t)=p(t)\rho(t)=p(t).