System utilities

  • order Order of a system.
  • size Number of outputs and inputs of a descriptor system .
  • iszero Checking whether the transfer function matrix of a descriptor system is zero.
  • evalfr Gain of the transfer function matrix at a single frequency value.
  • dcgain DC gain of a system.
  • opnorm L2- and L∞-norms of a descriptor system.
  • rss Generation of randomized standard state-space systems.
  • rdss Generation of randomized descriptor state-space systems.
  • dsxvarsel Building a descriptor systems by selecting a set of state variables.
  • dssubset Assigning a subsystem to a given descriptor system.
  • dszeros Setting a subsystem to zero.
  • dssubsel Selecting a subsystem according to a given zero-nonzero pattern.
  • dsdiag Building a k-times diagonal concatenation of a descriptor system.
DescriptorSystems.orderFunction
n = order(r)

Determine the order n of a rational transfer function r as the maximum of degrees of its numerator and denominator polynomials (n is also known as the McMillan degree of r).

source
nx = order(sys)

Return the order nx of the descriptor system sys as the dimension of the state variable vector. For improper or non-minimal systems, nx is less than the McMillan degree of the system.

source
Base.sizeFunction
size(sys) -> (p,m)
size(sys,1) -> p
size(sys,2) -> m

Return the number of outputs p and the number of inputs m of a descriptor system sys.

source
Base.iszeroFunction
 iszero(sys; atol = 0, atol1 = atol, atol2 = atol, rtol, fastrank = true)

Return true if the transfer function matrix of the descriptor system sys is zero. For a descriptor system sys = (A-λE,B,C,D) with the transfer function matrix G(λ) it is checked that the normal rank of G(λ) is zero, or equivalently (see [1]), that the normal rank of the system matrix pencil

          | A-λE | B | 
  S(λ) := |------|---|
          |  C   | D |

is equal to n, the order of the system sys.

If fastrank = true, the normal rank of S(λ) is evaluated by counting how many singular values of S(γ) have magnitudes greater than max(max(atol1,atol2), rtol*σ₁), where σ₁ is the largest singular value of S(γ) and γ is a randomly generated value. If fastrank = false, the rank is evaluated as nr + ni + nf + nl, where nr and nl are the sums of right and left Kronecker indices, respectively, while ni and nf are the number of infinite and finite eigenvalues, respectively. The sums nr+ni and nf+nl, are determined from an appropriate Kronecker-like form of the pencil S(λ), exhibiting the spliting of the right and left structures.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of A, B, C, D, the absolute tolerance for the nonzero elements of E, and the relative tolerance for the nonzero elements of A, B, C, D and E. The default relative tolerance is n*ϵ, where ϵ is the working machine epsilon. The keyword argument atol can be used to simultaneously set atol1 = atol and atol2 = atol.

[1] A. Varga, On checking null rank conditions of rational matrices, 2018. arXiv:2006.06825.

source
DescriptorSystems.evalfrFunction
Rval = evalfr(R,val)

Evaluate the rational transfer function matrix R(λ) for λ = val.

source
Rval = evalfr(R; fval = 0)

Evaluate the rational transfer function matrix R(λ) for λ = val, where val = im*fval for a continuous-time system or val = exp(im*fval*abs(Ts)) for a discrete-time system, with Ts the system sampling time.

source
rval = evalfr(r,val)

Evaluate the rational transfer function r(λ) for λ = val.

source
rval = evalfr(r; fval = 0)

Evaluate the rational transfer function r(λ) for λ = val, where val = im*fval for a continuous-time system or val = exp(im*fval*Ts) for a discrete-time system, with Ts the system sampling time.

source
Gval = evalfr(sys, val; atol1 = 0, atol2 = 0, rtol, fast = true)

Evaluate for a descriptor system sys = (A-λE,B,C,D) with the transfer function matrix G(λ), Gval, the value of the rational matrix G(λ) = C*inv(λE-A)*B+D for λ = val. The computed Gval has infinite entries if val is a pole (finite or infinite) of G(λ). If val is finite and val*E-A is singular or if val = Inf and E is singular, then the entries of Gval are evaluated separately for minimal realizations of each input-output channel.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices A, B, C, D, the absolute tolerance for the nonzero elements of E, and the relative tolerance for the nonzero elements of A, B, C, D and E.

The computation of minimal realizations of individual input-output channels relies on pencil manipulation algorithms, which employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

source
Gval = evalfr(sys; fval = 0, atol = 0, atol1 = atol, atol2 = atol, rtol, fast = true)

Evaluate for a descriptor system sys = (A-λE,B,C,D) with the transfer function matrix G(λ), Gval, the value of the rational matrix G(λ) = C*inv(λE-A)*B+D for λ = val, where val = im*fval for a continuous-time system or val = exp(im*fval*sys.Ts) for a discrete-time system. The computed Gval has infinite entries if val is a pole (finite or infinite) of G(λ). If val is finite and val*E-A is singular or if val = Inf and E is singular, then the entries of Gval are evaluated separately for minimal realizations of each input-output channel.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices A, B, C, D, the absolute tolerance for the nonzero elements of E, and the relative tolerance for the nonzero elements of A, B, C, D and E. The keyword argument atol can be used to simultaneously set atol1 = atol and atol2 = atol.

The computation of minimal realizations of individual input-output channels relies on pencil manipulation algorithms, which employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

source
DescriptorSystems.dcgainFunction
Gval = dcgain(sys; atol1, atol2, rtol, fast = true)

Evaluate for a descriptor system sys = (A-λE,B,C,D) with the transfer function matrix G(λ), Gval, the DC (or steady-state) gain. Gval is the value of the rational matrix G(λ) for λ = val, where for a continuous-time system val = 0 and for a discrete-time system val = 1. The computed Gval has infinite entries if val is a pole of G(λ). In this case (i.e., val*E-A is singular), the entries of Gval are evaluated separately for minimal realizations of each input-output channel.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices A, B, C, D, the absolute tolerance for the nonzero elements of E, and the relative tolerance for the nonzero elements of A, B, C, D and E.

The computation of minimal realizations of individual input-output channels relies on pencil manipulation algorithms, which employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

source
LinearAlgebra.opnormFunction
 opnorm(sys[, p = Inf]; kwargs...) 
 opnorm(sys, 2; kwargs...) -> sysnorm
 opnorm(sys, Inf; kwargs...) -> (sysnorm, fpeak)
 opnorm(sys; kwargs...) -> (sysnorm, fpeak)

Compute for a descriptor system sys = (A-λE,B,C,D) with the transfer function matrix G(λ), the L2 or L∞ system norm sysnorm induced by the vector p-norm, where valid values of p are 2 or Inf. For the L∞ norm, the frequency fpeak is also returned, where G(λ) achieves its peak gain. See gh2norm and ghinfnorm for a description of the allowed keyword arguments.

source
DescriptorSystems.rssFunction
sys = rss(n, p, m; disc = false, Ts, T = Float64, stable = false, nuc = 0, nuo = 0, randt = true)

Generate a randomized n+nuc+nuo-th order standard state-space system sys = (A,B,C,D) with p outputs and m inputs, with all matrices randomly generated of type T. The resulting sys is a continuous-time system if disc = false and a discrete-time system if disc = true. For a discrete-time system a sample time Δ can be specified using the keyword argument Ts = Δ (default: Δ = -1, i.e., not specified).

If stable = true, the resulting system is stable, with A having all eigenvalues with negative real parts for a continuous-time system, or with moduli less than one for a discrete-time system. If nuc+nuo > 0, the system sys is non-minimal, with A having nuc uncontrollable and nuo unobservable eigenvalues. If randt = true, a randomly generated orthogonal or unitary similarity transformation is additionally applied. If randt = false, the system matrices A, B, and C result in block stuctured forms exhibitting the uncontrollable and unobservable eigenvalues of A:

A = diag(A1, A2, A3),  B = [B1; 0; B3], C = [C1 C2 0]

with the diagonal blocks A1, A2, A3 of orders n, nuc, and nuo, respectively.

source
DescriptorSystems.rdssFunction
sys = rdss(n, p, m; id = [ ], disc = false, Ts, T = Float64, stable = false, nfuc = 0, iduc = [ ], 
           nfuo = 0, iduo = [ ], randlt = true, randrt = true)

Generate a randomized descriptor state-space system sys = (A-λE,B,C,D) with p outputs and m inputs, with all matrices randomly generated of type T. The resulting sys is a continuous-time system if disc = false and a discrete-time system if disc = true. For a discrete-time system a sample time Δ can be specified using the keyword argument Ts = Δ (default: Δ = -1, i.e., not specified).

If the vector id is nonempty, then id[i] specifies the order of the i-th infinite elementary divisor of the resulting pencil A-λE, which thus has n finite eigenvalues and ni = sum(id) infinite eigenvalues which are controllable and observable. If nfuc+nfuo > 0, the system sys is non-minimal, with A having nfuc uncontrollable and nfuo unobservable finite eigenvalues. If the vector iduc is a nonempty, then iduc[i] specifies the order of the i-th infinite elementary divisor with uncontrollable infinite eigenvalues of the resulting pencil A-λE, which thus has niuc = sum(iduc) uncontrollable infinite eigenvalues. If the vector iduo is a nonempty, then iduo[i] specifies the order of the i-th infinite elementary divisor with unobservable infinite eigenvalues of the resulting pencil A-λE, which thus has niuo = sum(iduo) unobservable infinite eigenvalues. If niuc+niuo > 0, the system sys is non-minimal, with A having niuc uncontrollable and niuo unobservable infinite eigenvalues.

It follows, that the resulting pencil A-λE has n+nfuc+nfuo finite eigenvalues and ni+niuc+niuo infinite eigenvalues. If stable = true, the proper part of the system sys is stable, with A having all finite eigenvalues with negative real parts for a continuous-time system, or with moduli less than one for a discrete-time system.

If randlt = true, a randomly generated orthogonal or unitary transformation is additionally applied to A, E, and B from the left. If randrt = true, a randomly generated orthogonal or unitary transformation is additionally applied to A, E, and C from the right. If randlt = false and randrt = false, the system matrices A, E, B, and C result in block stuctured forms exhibitting the uncontrollable and unobservable finite and infinite eigenvalues of A-λE:

A-λE = diag(A1-λE1, A2-λE2, A3-λE3, A4-λE4, A5-λE5, A6-λE6),  
B = [B1; B2; 0; 0; B5; B6 ], 
C = [C1 C2 C3 C4 0 0]

with the diagonal blocks A1, A2, A3, A4, A5, A6 of orders n, ni, nfuc, niuc, nfuo and niuo, respectively.

source
DescriptorSystems.dsxvarselFunction
sysr = dsxvarsel(sys,ind)

Construct for the descriptor system sys = (A-λE,B,C,D) of order n the descriptor system sysr = (A[ind,ind]-λE[ind,ind],B[ind,:],C[:,ind],D) of order nr = length(ind), by selecting the state variables of sys with indices specified by ind. If ind is a permutation vector of length n, then sysr has the same transfer function matrix as sys and permuted state variables.

source
DescriptorSystems.dssubsetFunction
 res = dssubset(sys, subsys, rows, cols; minimal = true, fast = true, atol = 0, atol1 = atol, atol2 = atol, rtol = nϵ)

Assign for a given descriptor system sys = (A-λE,B,C,D) its subsystem sys[row,cols] to subsys and return the modified system in res. rows and cols are indices, vectors of indices, index ranges, : or any combinations of them.

If minimal = true (default), an irreducible realization of the resulting system res is computed, otherwise a possibly non-minimal realization is returned if minimal = false.

The underlying pencil manipulation based minimal realization algorithms employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices A, B, C, D, the absolute tolerance for the nonzero elements of E, and the relative tolerance for the nonzero elements of A, B, C, D and E. The default relative tolerance is , where ϵ is the working machine epsilon and n is the order of the system sys. The keyword argument atol can be used to simultaneously set atol1 = atol and atol2 = atol.

source
DescriptorSystems.dszerosFunction
 res = dszeros(sys, rows, cols; minimal = true, fast = true, atol = 0, atol1 = atol, atol2 = atol, rtol = nϵ)

Set for a given descriptor system sys = (A-λE,B,C,D) its subsystem sys[row,cols] to zero and return the modified system in res. rows and cols are indices, vectors of indices, index ranges, : or any combinations of them.

If minimal = true (default), an irreducible realization of the resulting system res is computed, otherwise a possibly non-minimal realization is returned if minimal = false.

The underlying pencil manipulation based minimal realization algorithms employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices A, B, C, D, the absolute tolerance for the nonzero elements of E, and the relative tolerance for the nonzero elements of A, B, C, D and E. The default relative tolerance is , where ϵ is the working machine epsilon and n is the order of the system sys. The keyword argument atol can be used to simultaneously set atol1 = atol and atol2 = atol.

source
DescriptorSystems.dssubselFunction
 res = dssubsel(sys, S; minimal = true, fast = true, atol = 0, atol1 = atol, atol2 = atol, rtol = nϵ)

Select for a given descriptor system sys = (A-λE,B,C,D) its subsystem res corresponding to the zero-one pattern specified by the binary structure matrix S. If G(λ) is the transfer function matrix of sys, then the transfer function matrix of res is the element-wise product S .* G(λ).

If minimal = true (default), an irreducible realization of the resulting system res is computed, otherwise a possibly non-minimal realization is returned if minimal = false.

The underlying pencil manipulation based minimal realization algorithms employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true, or the SVD-decomposition. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of matrices A, B, C, D, the absolute tolerance for the nonzero elements of E, and the relative tolerance for the nonzero elements of A, B, C, D and E. The default relative tolerance is , where ϵ is the working machine epsilon and n is the order of the system sys. The keyword argument atol can be used to simultaneously set atol1 = atol and atol2 = atol.

source
DescriptorSystems.dsdiagFunction
 sysdiag = dsdiag(sys, k)

Build for a given descriptor system sys and non-negative integer k a descriptor system sysdiag obtained as sysdiag = diag( sys, ..., sys ) (i.e., k-times repeated application of append).

source