Building descriptor system state-space models

  • DescriptorStateSpace Descriptor state-space object.
  • dss Construction of descriptor state-space models.
  • dssdata Extraction of matrix data from a descriptor state-space model.
DescriptorSystems.DescriptorStateSpaceType
DescriptorStateSpace{T}(A::Matrix{T}, E::Union{Matrix{T},UniformScaling}, 
                        B::Matrix{T}, C::Matrix{T}, D::Matrix{T},  
                        Ts::Real) where T <: Number

Construct a descriptor state-space model from a quintuple of matrices (A,E,B,C,D) and a sampling time Ts.

If SYS::DescriptorStateSpace{T} is a descriptor system model object defined by the 4-tuple SYS = (A-λE,B,C,D), then:

SYS.A is the nx × nx state matrix A with elements of type T.

SYS.E is the nx × nx descriptor matrix E with elements of type T. For a standard state-space system SYS.E = I, the UniformScaling of type Bool.

SYS.B is the nx × nu system input matrix B with elements of type T.

SYS.C is the ny × nx system output matrix C with elements of type T.

SYS.D is the ny × nu system feedthrough matrix D with elements of type T.

SYS.Ts is the real sampling time Ts, where Ts = 0 for a continuous-time system, and Ts > 0 or Ts = -1 for a discrete-time system. Ts = -1 indicates a discrete-time system with an unspecified sampling time.

The dimensions nx, ny and nu can be obtained as SYS.nx, SYS.ny and SYS.nu, respectively.

source
DescriptorSystems.dssFunction
sys = dss(A, E, B, C, D; Ts = 0, check_reg = false, 
          atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ )

Create for Ts = 0 a descriptor system model sys = (A-λE,B,C,D) for a continuous-time state space system of the form

Edx(t)/dt = Ax(t) + Bu(t) ,
y(t)      = Cx(t) + Du(t) ,

where x(t), u(t) and y(t) are the system state vector, system input vector and system output vector, respectively, for the continuous time variable t.

For a nonzero positive sampling time Ts = ΔT, the descriptor system model specifies a discrete-time state space system of the form

Ex(t+ΔT) = Ax(t) + Bu(t)
y(t)     = Cx(t) + Du(t)

for the discrete values of the time variable t = 0, ΔT, 2ΔT, .... Use Ts = -1 if the sampling time is not specified. In this case, by convention ΔT = 1.

For a system with zero feedthrough matrix D, it is possible to set D = 0 (the scalar zero).

For a standard state space system, E is the identity matrix. In this case, it is possible to set E = I (the boolean uniform scaling). Alternatively, use

sys = dss(A, B, C, D; Ts = 0)

to create a standard system.

For a system corresponding to a static gain D, use

sys = dss(D; Ts = 0)

It is possible to specify a descriptor system via all or part of its matrices using the form

sys = dss(A = mat1, E = mat2, B = mat3, C = mat4, D = mat5; Ts = 0, check_reg = false, 
          atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ)

where A, E, B, C, and D are keyword parameters set to appropriate matrix values mat1, mat2, mat3, mat4, and mat5, respectively. If some of the system matrices are omited, then zero matrices of appropriate sizes are employed instead.

It is assumed that the pencil A-λE is regular (i.e., det(A-λE) ̸≡ 0), and therefore, in the interest of efficiency, the regularity of A-λE is by default not tested. If check_reg = true, the regularity of A-λE is additionally checked. In this case, the keyword arguments atol1, atol2 and rtol specify, respectively, the absolute tolerance for the nonzero elements of A, the absolute tolerance for the nonzero elements of E, and the relative tolerance for the nonzero elements of A and E. The default relative tolerance is n*ϵ, where n is the order of the square matrices A and E, and ϵ is the working machine epsilon. The keyword argument atol can be used to simultaneously set atol1 = atol and atol2 = atol.

source
 sys = dss(A, E, B, F, C, G, D, H; compacted = false, 
           atol = 0, atol1 = atol, atol2 = atol, atol3 = atol, rtol = min(atol1,atol2,atol3)>0 ? 0 : n*ϵ)

Construct an input-output equivalent descriptor system representation sys = (Ad-λdE,Bd,Cd,Dd) to a pencil based linearization (A-λE,B-λF,C-λG,D-λH) satisfying

            -1                        -1
 Cd*(λEd-Ad)  *Bd + Dd = (C-λG)*(λE-A)  *(B-λF) + D-λH .

If compacted = true, a compacted descriptor system realization is determined by exploiting possible rank defficiencies of the matrices F, G, and H. Any of the matrices F, G, and H can be set to missing.

The keyword arguments atol1, atol2, atol3, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of F, the absolute tolerance for the nonzero elements of G, the absolute tolerance for the nonzero elements of H and the relative tolerance for the nonzero elements of F, G and H. The default relative tolerance is n*ϵ, where n is the size of of A, and ϵ is the machine epsilon of the element type of A. The keyword argument atol can be used to simultaneously set atol1 = atol, atol2 = atol and atol3 = atol.

source
sys = dss(NUM, DEN; contr = false, obs = false, noseig = false, minimal = false, fast = true, atol = 0, rtol)

Convert the rational matrix R(λ) = NUM(λ) ./ DEN(λ) to a descriptor system representation sys = (A-λE,B,C,D) such that the transfer function matrix of sys is R(λ).

NUM(λ) is a polynomial matrix of the form NUM(λ) = N_1 + λ N_2 + ... + λ**k N_(k+1), for which the coefficient matrices N_i, i = 1, ..., k+1 are stored in the 3-dimensional matrix NUM, where NUM[:,:,i] contains the i-th coefficient matrix N_i (multiplying λ**(i-1)).

DEN(λ) is a polynomial matrix of the form DEN(λ) = D_1 + λ D_2 + ... + λ**l D_(l+1), for which the coefficient matrices D_i, i = 1, ..., l+1, are stored in the 3-dimensional matrix DEN, where DEN[:,:,i] contain the i-th coefficient matrix D_i (multiplying λ**(i-1)).

Alternatively, NUM(λ) and DEN(λ) can be specified as matrices of elements of the Polynomial type provided by the Polynomials package.

If n is the order of A-λE, then the computed linearization satisfies:

(1) A-λE is regular and R(λ) = C*inv(λE-A)*B+D;

(2) rank[B A-λE] = n (controllability) if minimal = true or contr = true;

(3) rank[A-λE; C] = n (observability) if minimal = true or obs = true;

(4) A-λE has no non-dynamic modes if minimal = true or noseig = true.

If conditions (1)-(4) are satisfied, the realization is called minimal and the resulting order n is the least achievable order. If conditions (1)-(3) are satisfied, the realization is called irreducible and the resulting order n is the least achievable order using orthogonal similarity transformations. An irreducible realization preserves the pole-zero and singular structures of R(λ).

The descriptor system based realization is built using the methods described in [1] in conjunction with pencil manipulation algorithms [2] and [3] to compute reduced order realization. These algorithms employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition, if fast = false. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

The keyword arguments atol and rtol, specify the absolute and relative tolerances, respectively, for the nonzero coefficients of NUM(λ) and DEN(λ).

[1] A. Varga, On computing the Kronecker structure of polynomial and rational matrices using Julia, 2020, arXiv:2006.06825.

[2] P. Van Dooreen, The generalized eigenstructure problem in linear system theory, IEEE Transactions on Automatic Control, vol. AC-26, pp. 111-129, 1981.

[3] A. Varga, Solving Fault Diagnosis Problems - Linear Synthesis Techniques, Springer Verlag, 2017.

source
sys = dss(R; Ts=missing, contr = false, obs = false, noseig = false, minimal = false, fast = true, atol = 0, rtol)

Convert the rational transfer function matrix R(λ) to a descriptor system representation sys = (A-λE,B,C,D) such that the transfer function matrix of sys is R(λ). The resulting sys is a continuous-time system if Ts = 0 or discrete-time system if Ts = -1 or Ts > 0. If Ts = missing, the sampling time of sys is inherited from the sampling time TRs of the elements of R, unless TRs = nothing, in which case Ts = 0 is used (by default).

R(λ) is a matrix with rational transfer function entries (see RationalTransferFunction ) corresponding to a multiple-input multiple-outputs system or a rational transfer function corresponding to a single-input single-output system. The numerators and denominators of the elements of R are of type Polynomial as provided by the Polynomials package.

If n is the order of A-λE, then the computed realization satisfies:

(1) A-λE is regular and R(λ) = C*inv(λE-A)*B+D;

(2) rank[B A-λE] = n (controllability) if minimal = true or contr = true;

(3) rank[A-λE; C] = n (observability) if minimal = true or contr = true;

(4) A-λE has no non-dynamic modes if minimal = true or noseig = true.

If conditions (1)-(4) are satisfied, the realization is called minimal and the resulting order n is the least achievable order. If conditions (1)-(3) are satisfied, the realization is called irreducible and the resulting order n is the least achievable order using orthogonal similarity transformations. An irreducible realization preserves the pole-zero and singular structures of R(λ).

The underlying pencil manipulation algorithms [1] and [2] to compute reduced order realizations employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition, if fast = false. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

The keyword arguments atol and rtol, specify the absolute and relative tolerances for the nonzero coefficients of R(λ).

[1] P. Van Dooreen, The generalized eigenstructure problem in linear system theory, IEEE Transactions on Automatic Control, vol. AC-26, pp. 111-129, 1981.

[2] A. Varga, Solving Fault Diagnosis Problems - Linear Synthesis Techniques, Springer Verlag, 2017.

source
sys = dss(P; Ts = 0, contr = false, obs = false, noseig = false, minimal = false, fast = true, atol = 0, rtol)

Convert the polynomial matrix P(λ) to a descriptor system representation sys = (A-λE,B,C,D) such that the transfer function matrix of sys is P(λ). The resulting sys is a continuous-time system if Ts = 0 or discrete-time system if Ts = -1 or Ts > 0.

P(λ) can be specified as a grade k polynomial matrix of the form P(λ) = P_1 + λ P_2 + ... + λ**k P_(k+1), for which the coefficient matrices P_i, i = 1, ..., k+1, are stored in the 3-dimensional matrix P, where P[:,:,i] contains the i-th coefficient matrix P_i (multiplying λ**(i-1)).

P(λ) can also be specified as a matrix, vector or scalar of elements of the Polynomial type provided by the Polynomials package.

If d is the degree of P(λ) and n is the order of A-λE, then the computed realization satisfies:

(1) A-λE is regular and P(λ) = C*inv(λE-A)*B+D;

(2) rank[B A-λE] = n (controllability) if minimal = true or contr = true;

(3) rank[A-λE; C] = n (observability) if minimal = true or contr = true;

(4) A-λE has no non-dynamic modes if minimal = true or noseig = true.

If conditions (1)-(4) are satisfied, the realization is called minimal and the resulting order n is the least achievable order. If conditions (1)-(3) are satisfied, the realization is called irreducible and the resulting order n is the least achievable order using orthogonal similarity transformations. An irreducible realization preserves the pole-zero and singular structures of P(λ).

The underlying pencil manipulation algorithms [1] and [2] to compute reduced order realizations employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition, if fast = false. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

The keyword arguments atol and rtol, specify the absolute and relative tolerances for the nonzero coefficients of P(λ), respectively.

[1] P. Van Dooreen, The generalized eigenstructure problem in linear system theory, IEEE Transactions on Automatic Control, vol. AC-26, pp. 111-129, 1981.

[2] A. Varga, Solving Fault Diagnosis Problems - Linear Synthesis Techniques, Springer Verlag, 2017.

source
sys = dss(T, U, V, W; fast = true, contr = false, obs = false, minimal = false, atol = 0, rtol)

Construct an input-output equivalent descriptor system representation sys = (A-λE,B,C,D) to a polynomial model specified by the polynomial matrices T(λ), U(λ), V(λ), and W(λ) such that

  V(λ)*inv(T(λ))*U(λ)+W(λ) = C*inv(λE-A)*B+D.

If minimal = true, the resulting realization (A-λE,B,C,D) has the least possible order n of A-λE.

T(λ), U(λ), V(λ), and W(λ) can be specified as polynomial matrices of the form X(λ) = X_1 + λ X_2 + ... + λ**k X_(k+1), for X = T, U, V, and W, for which the coefficient matrices X_i, i = 1, ..., k+1, are stored in the 3-dimensional matrices X, where X[:,:,i] contains the i-th coefficient matrix X_i (multiplying λ**(i-1)).

T(λ), U(λ), V(λ), and W(λ) can also be specified as matrices, vectors or scalars of elements of the Polynomial type provided by the Polynomials package. In this case, no check is performed that T(λ), U(λ), V(λ) and W(λ) have the same indeterminates.

The computed descriptor realization satisfies:

(1) A-λE is regular;

(2) rank[B A-λE] = n (controllability) if minimal = true or contr = true;

(3) rank[A-λE; C] = n (observability) if minimal = true or obs = true;

(4) A-λE has no simple infinite eigenvalues if minimal = true.

The keyword arguments atol and rtol, specify, respectively, the absolute and relative tolerance for the nonzero coefficients of the matrices T(λ), U(λ), V(λ) and W(λ). The default relative tolerance is nt*ϵ, where nt is the size of the square matrix T(λ) and ϵ is the machine epsilon of the element type of its coefficients.

The descriptor realization is built using the methods described in [1].

[1] A. Varga, On computing the Kronecker structure of polynomial and rational matrices using Julia, 2020, arXiv:2006.06825.

source
DescriptorSystems.dssdataFunction
A, E, B, C, D  = dssdata([T,] sys)

Extract the matrices A, E, B, C, D of a descriptor system model sys = (A-λE,B,C,D). If the type T is specified, the resulting matrices are converted to this type.

source