Descriptor system analysis
isregular
Test whether a descriptor system has a regular pole pencil.gpole
Poles of a descriptor system.gpoleinfo
Poles and pole structure information of a descriptor system.isproper
Test whether a descriptor system is proper.isstable
Test whether a descriptor system is stable.gzero
Zeros of a descriptor system.gzeroinfo
Zeros and zero structure information of a descriptor system.gnrank
Normal rank of the transfer function matrix of a descriptor system.ghanorm
Hankel norm of a proper and stable descriptor system.gl2norm
L2
norm of a descriptor system.gh2norm
H2
norm of a descriptor system.glinfnorm
L∞
norm of a descriptor system.ghinfnorm
H∞
norm of a descriptor system.gnugap
ν-gap
distance between two descriptor systems.freqresp
Frequency response of a descriptor system.timeresp
Time response of a descriptor system.stepresp
Step response of a descriptor system.gbalqual
Evaluation of the scaling quality of the matrices of a descriptor system.
MatrixPencils.isregular
— Functionisregular(sys; atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ)
Return true
if the descriptor system sys = (A-λE,B,C,D)
has a regular pole pencil A-λE
and false
otherwise.
To test whether the pencil A-λE
is regular (i.e., det(A-λE) ̸≡ 0
), the underlying computational procedure reduces the pencil A-λE
to an appropriate Kronecker-like form, which provides information on the rank of A-λE
.
The keyword arguements atol1
, atol2
and rtol
specify 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
, respectively. The default relative tolerance is n*ϵ
, where n
is the size of A
, and ϵ
is the working machine epsilon. The keyword argument atol
can be used to simultaneously set atol1 = atol
and atol2 = atol
.
DescriptorSystems.gpole
— Functionval = gpole(sys; fast = false, atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ)
Return for the descriptor system sys = (A-λE,B,C,D)
the complex vector val
containing the finite and infinite zeros of the system pole pencil P(λ) := A-λE
. The values in val
are the poles of the transfer function matrix of sys
, if A-λE
is regular and the descriptor system realization sys = (A-λE,B,C,D)
is irreducible. If the pencil A-λE
is singular, val
also contains NaN
elements, whose number is the rank deficiency of the pencil A-λE
.
For E
nonsingular, val
contains the generalized eigenvalues of the pair (A,E)
. For E
singular, val
contains the zeros of P(λ)
, which are computed by reducing the pencil P(λ)
to an appropriate Kronecker-like form using orthonal similarity transformations and involves rank decisions based on rank revealing QR-decompositions with column pivoting, if fast = true
, or, the more reliable, SVD-decompositions, if fast = false
.
The regularity of A-λE
is implicitly checked. If check_reg = true
, an error message is issued if the pencil A-λE
is singular. If check_reg = false
and the pencil A-λE
is singular, then n-r
poles are set to NaN
, where n
is the system order and r
is the normal rank of A-λE
.
The keyword arguements atol1
, atol2
and rtol
specify 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
, respectively. 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
.
DescriptorSystems.gpoleinfo
— Functiongpoleinfo(sys; smarg, fast = false, atol = 0, atol1 = atol, atol2 = atol,
rtol = n*ϵ, offset = sqrt(ϵ)) -> (val, info)
Return for the descriptor system sys = (A-λE,B,C,D)
the complex vector val
containing the finite and infinite zeros of the system pole pencil P(λ) := A-λE
and the named tuple info
containing information on the eigenvalue structure of the pole pencil P(λ)
. The values in val
are the poles of the transfer function matrix of sys
, if A-λE
is regular and the descriptor system realization sys = (A-λE,B,C,D)
is irreducible. If the pencil A-λE
is singular, val
also contains NaN
elements, whose number is the rank deficiency of the pencil A-λE
.
For stability analysis purposes, a stability margin smarg
can be specified for the finite eigenvalues, in conjunction with a stability domain boundary offset β
to numerically assess the finite eigenvalues which belong to the boundary of the stability domain as follows: in the continuous-time case, these are the finite eigenvalues having real parts in the interval [smarg-β, smarg+β]
, while in the discrete-time case, these are the finite eigenvalues having moduli in the interval [smarg-β, smarg+β]
. The default value of the stability margin smarg
is 0
for a continuous-time system and 1
for a discrete-time system. The default value used for β
is sqrt(ϵ)
, where ϵ
is the working machine precision.
The named tuple info
contains the following information:
info.nfev
is the number of finite eigenvalues of the pencil A-λE
(also the number of finite poles of sys
);
info.niev
is the number of infinite eigenvalues of the pencil A-λE
;
info.nisev
is the number of simple infinite eigenvalues of the pencil A-λE
(also known as non-dynamic modes);
info.nip
is the number of infinite poles of the system sys
;
info.nfsev
is the number of finite stable eigenvalues, i.e., the finite eigenvalues having real parts or moduli less than smarg-β
for a continuous- or discrete-time system, respectively;
info.nfsbev
is the number of finite eigenvalues on the boundary of the stability domain, i.e., the finite eigenvalues having real parts or moduli in the interval [smarg-β, smarg+β]
for a continuous- or discrete-time system, respectively;
info.nfuev
is the number of finite unstable eigenvalues, i.e., the finite eigenvalues having real parts or moduli greater than smarg+β
for a continuous- or discrete-time system, respectively;
info.nhev
is the number of hidden eigenvalues set to NaN
(can be nonzero only if the pencil A-λE
is singular);
info.nrank
is the normal rank of the pencil A-λE
;
info.miev
is an integer vector, which contains the multiplicities of the infinite eigenvalues of the pencil A-λE
as follows: the i
-th element info.miev[i]
is the order of an infinite elementary divisor (i.e., the multiplicity of an infinite eigenvalue) and the number of infinite poles is the sum of the components of info.miev
;
info.mip
is an integer vector, which contains the information on the multiplicities of the infinite zeros of A-λE
as follows: the i
-th element info.mip[i]
is equal to k-1
, where k
is the order of an infinite elementary divisor with k > 0
and the number of infinite poles is the sum of the components of info.mip
;
info.rki
is an integer vector, which contains the right Kronecker indices of the pencil A-λE
(empty for a regular pencil);
info.lki
is an integer vector, which contains the left Kronecker indices of the pencil A-λE
(empty for a regular pencil);
info.regular
is set to true
, if the pencil A-λE
is regular and set to false
, if the pencil A-λE
is singular;
info.proper
is set to true
, if the pencil A-λE
is regular and all its infinite eigenvalues are simple (has only non-dynamic modes), or is set to false
, if the pencil A-λE
is singular or has higher order infinite eigenvalues;
info.stable
is set to true
, if the pencil A-λE
is regular, has only stable finite eigenvalues and all its infinite eigenvalues are simple (has only non-dynamic modes), and is set to false
otherwise.
Note: The finite poles and the finite eigenvalues of the pencil P(λ)
are the same, but the multiplicities of infinite eigenvalues of P(λ)
are in excess with one to the multiplicities of infinite poles.
For the reduction of the pencil P(λ)
to an appropriate Kronecker-like form orthonal similarity transformations are performed, which involve rank decisions based on rank revealing QR-decompositions with column pivoting, if fast = true
, or, the more reliable, SVD-decompositions, if fast = false
.
The keyword arguements atol1
, atol2
and rtol
specify 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
, respectively. The default relative tolerance is n*ϵ
, where n
is the size of P(λ)
, and ϵ
is the working machine epsilon. The keyword argument atol
can be used to simultaneously set atol1 = atol
and atol2 = atol
.
DescriptorSystems.isproper
— Functionisproper(sys; atol = 0, atol1 = atol, atol2 = atol, rtol = = n*ϵ, fast = true)
Return true
if the transfer function matrix G(λ)
of the descriptor system sys = (A-λE,B,C,D)
is proper and false
otherwise.
For a descriptor system realization sys = (A-λE,B,C,D)
without uncontrollable and unobservable infinite eigenvalues, it is checked that the pencil A-λE
has no infinite eigenvalues or, if infinite eigenvalues exist, all infinite eigenvalues are simple. If the original descriptor realization has uncontrollable or unobservable infinite eigenvalues, these are elliminated using orthogonal pencil reduction algorithms.
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 n
is the order of A
and ϵ
is the working machine epsilon. The keyword argument atol
can be used to simultaneously set atol1 = atol
and atol2 = atol
.
DescriptorSystems.isstable
— Functionisstable(sys[, smarg]; fast = true, atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ, offset = sqrt(ϵ))
Return true
if the descriptor system sys = (A-λE,B,C,D)
has only stable poles and false
otherwise.
It is checked that the pole pencil P(λ) := A-λE
has no infinite eigenvalues or, if infinite eigenvalues exist, all infinite eigenvalues are simple, and additionally the real parts of all finite eigenvalues are less than smarg-β
for a continuous-time system or have moduli less than smarg-β
for a discrete-time system, where smarg
is the stability margin and β
is the stability domain boundary offset. The default value of the stability margin smarg
is 0
for a continuous-time system and 1
for a discrete-time system. The offset β
to be used to numerically assess the stability of eigenvalues can be specified via the keyword parameter offset = β
. The default value used for β
is sqrt(ϵ)
, where ϵ
is the working machine precision.
For E
singular, the computation of the poles is performed by reducing the pencil P(λ)
to an appropriate Kronecker-like form using orthonal similarity transformations and involves rank decisions based on rank revealing QR-decompositions with column pivoting, if fast = true
, or, the more reliable, SVD-decompositions, if fast = false
.
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 n
is the order of A
and ϵ
is the working machine epsilon. The keyword argument atol
can be used to simultaneously set atol1 = atol
and atol2 = atol
.
DescriptorSystems.gzero
— Functionval = gzero(sys; fast = false, prescale, atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ)
Return for the descriptor system sys = (A-λE,B,C,D)
the complex vector val
containing the finite and infinite Smith zeros of the system matrix pencil
| A-λE | B |
S(λ) := |------|---| .
| C | D |
The values in val
are called the invariant zeros of the pencil S(λ)
and are the transmission zeros of the transfer function matrix of sys
if A-λE
is regular and the descriptor system realization sys = (A-λE,B,C,D)
is irreducible.
If prescale = true
, a preliminary balancing of the descriptor system matrices is performed. The default setting is prescale = gbalqual(sys) > 10000
, where gbalqual(sys)
is the scaling quality of the descriptor system model sys
(see gbalqual
).
The computation of the zeros is performed by reducing the pencil S(λ)
to an appropriate Kronecker-like form using orthonal similarity transformations and involves rank decisions based on rank revealing QR-decompositions with column pivoting, if fast = true
, or, the more reliable, SVD-decompositions, if fast = false
.
The keyword arguements atol1
, atol2
and rtol
specify the absolute tolerance for the nonzero elements of A
, B
, C
and D
, the absolute tolerance for the nonzero elements of E
, and the relative tolerance for the nonzero elements of A
, E
, B
, C
and D
, respectively. The default relative tolerance is n*ϵ
, where n
is the size of A
, and ϵ
is the working machine epsilon. The keyword argument atol
can be used to simultaneously set atol1 = atol
and atol2 = atol
.
DescriptorSystems.gzeroinfo
— Functiongzeroinfo(sys; smarg, fast = false, prescale, atol = 0, atol1 = atol, atol2 = atol,
rtol = n*ϵ, offset = sqrt(ϵ)) -> (val, info)
Return for the descriptor system sys = (A-λE,B,C,D)
the complex vector val
containing the finite and infinite Smith zeros of the system matrix pencil S(λ)
| A-λE | B |
S(λ) = |------|---|
| C | D |
and the named tuple info
containing information on the Kronecker structure of the pencil S(λ)
. The values in val
are called the invariant zeros of the pencil S(λ)
and are the transmission zeros of the transfer function matrix of sys
if A-λE
is regular and the descriptor system realization sys = (A-λE,B,C,D)
is irreducible.
If prescale = true
, a preliminary balancing of the descriptor system matrices is performed. The default setting is prescale = gbalqual(sys) > 10000
, where gbalqual(sys)
is the scaling quality of the descriptor system model sys
(see gbalqual
).
For stability analysis purposes, a stability margin smarg
can be specified for the finite zeros, in conjunction with a stability domain boundary offset β
to numerically assess the finite zeros which belong to the boundary of the stability domain as follows: in the continuous-time case, these are the finite zeros having real parts in the interval [smarg-β, smarg+β]
, while in the discrete-time case, these are the finite zeros having moduli in the interva [smarg-β, smarg+β]
. The default value of the stability margin smarg
is 0
for a continuous-time system and 1
for a discrete-time system. The default value used for β
is sqrt(ϵ)
, where ϵ
is the working machine precision.
The named tuple info
contains the following information:
info.nfz
is the number of finite eigenvalues of the pencil S(λ)
(also the number of finite zeros of sys
);
info.niev
is the number of infinite eigenvalues of the pencil S(λ)
;
info.nisev
is the number of simple infinite eigenvalues of the pencil S(λ)
;
info.niz
is the number of infinite zeros of the system sys
;
info.nfsz
is the number of finite stable zeros, i.e., the finite zeros having real parts or moduli less than smarg-β
for a continuous- or discrete-time system, respectively;
info.nfsbz
is the number of finite zeros on the boundary of the stability domain, i.e., the finite zeros having real parts or moduli in the interval [smarg-β, smarg+β]
for a continuous- or discrete-time system, respectively;
info.nfuz
is the number of finite unstable zeros, i.e., the finite zeros having real parts or moduli greater than smarg+β
for a continuous- or discrete-time system, respectively;
info.nrank
is the normal rank of the pencil S(λ)
;
info.miev
is an integer vector, which contains the multiplicities of the infinite eigenvalues of the pencil S(λ)
(also the dimensions of the elementary infinite blocks in the Kronecker form of S(λ)
);
info.miz
is an integer vector, which contains the information on the multiplicities of the infinite zeros of S(λ)
as follows: S(λ)
has info.mip[i]
infinite zeros of multiplicity i
, and is empty if S(λ)
has no infinite zeros;
info.rki
is an integer vector, which contains the right Kronecker indices of the pencil S(λ)
(empty for a regular pencil);
info.lki
is an integer vector, which contains the left Kronecker indices of the pencil S(λ)
(empty for a regular pencil);
info.regular
is set to true
, if the pencil S(λ)
is regular and set to false
, if the pencil S(λ)
is singular;
info.stable
is set to true
, if the pencil S(λ)
has only stable finite zeros and all its infinite zeros are simple and is set to false
otherwise.
Note: The finite zeros and the finite eigenvalues of the pencil S(λ)
are the same, but the multiplicities of infinite eigenvalues are in excess with one to the multiplicities of infinite zeros.
The computation of the zeros is performed by reducing the pencil S(λ)
to an appropriate Kronecker-like form using orthonal similarity transformations and involves rank decisions based on rank revealing QR-decompositions with column pivoting, if fast = true
, or, the more reliable, SVD-decompositions, if fast = false
.
The keyword arguements atol1
, atol2
and rtol
specify the absolute tolerance for the nonzero elements of A
, B
, C
and D
, the absolute tolerance for the nonzero elements of E
, and the relative tolerance for the nonzero elements of A
, E
, B
, C
and D
, respectively. The default relative tolerance is n*ϵ
, where n
is the size of A
and ϵ
is the working machine epsilon. The keyword argument atol
can be used to simultaneously set atol1 = atol
and atol2 = atol
.
DescriptorSystems.gnrank
— Functionr = gnrank(sys, fastrank = true, atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ )
Compute the normal rank r
of the transfer function matrix G(λ)
of the descriptor system sys = (A-λE,B,C,D)
.
The normal rank of G(λ)
is evaluated as r = k - n
, where k
is the normal rank of the system matrix pencil
| A-λE | B |
S(λ) := |------|---|
| C | D |
and n
is the order of the system sys
(i.e., the size of A
).
If fastrank = true
, the normal rank of S(λ)
is evaluated by counting the singular values of S(γ)
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
.
DescriptorSystems.ghanorm
— Functionghanorm(sys, fast = true, atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ) -> (hanorm, hs)
Compute for a proper and stable descriptor system sys = (A-λE,B,C,D)
with the transfer function matrix G(λ)
, the Hankel norm hanorm =
$\small ||G(\lambda)||_H$ and the vector of Hankel singular values hs
of the minimal realizatioj of the system.
For a non-minimal system, the uncontrollable and unobservable finite and infinite eigenvalues of the pair (A,E)
and the non-dynamic modes are elliminated using minimal realization techniques. The rank determinations in the performed reductions are based on rank revealing QR-decompositions with column pivoting if fast = true
or the more reliable SVD-decompositions if fast = false
.
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 and n
is the order of the system sys
. The keyword argument atol
can be used to simultaneously set atol1 = atol
and atol2 = atol
.
DescriptorSystems.gl2norm
— Functiongl2norm(sys, h2norm = false, fast = true, offset = sqrt(ϵ), atol = 0, atol1 = atol, atol2 = atol, atol3 = atol, atolinf = atol, rtol = n*ϵ)
Compute for a descriptor system sys = (A-λE,B,C,D)
the L2
norm of its transfer function matrix G(λ)
. The L2
norm is infinite if G(λ)
has poles on the stability domain boundary, i.e., on the extended imaginary axis, in the continuous-time case, or on the unit circle, in the discrete-time case. The L2
norm is also infinite for a continuous-time system having a gain at infinity greater than atolinf
. If the pencil A-λE
has uncontrollable and/or unobservable eigenvalues on the boundary of the stability domain, then a reduced order realization is determined first (see below) to eliminate these eigenvalues.
To check the lack of poles on the stability domain boundary, the eigenvalues of the pencil A-λE
must not have real parts in the interval [-β,β]
for a continuous-time system or must not have moduli in the interval [1-β,1+β]
for a discrete-time system, where β
is the stability domain boundary offset. The offset β
to be used can be specified via the keyword parameter offset = β
. The default value used for β
is sqrt(ϵ)
, where ϵ
is the working machine precision.
If h2norm = true
, the H2
norm is computed. The H2
norm is infinite if the pole pencil A-λE
has unstable zeros (i.e., unstable poles), or for a continuous-time system having a gain at infinity greater than atolinf
. To check the stability, the eigenvalues of the pencil A-λE
must have real parts less than -β
for a continuous-time system or have moduli less than 1-β
for a discrete-time system.
For a continuous-time system sys
with E
singular, a reduced order realization is determined first, without uncontrollable and unobservable finite and infinite eigenvalues of the pencil A-λE
. For a discrete-time system or for a system with invertible E
, a reduced order realization is determined first, without uncontrollable and unobservable finite eigenvalues of the pencil A-λE
. The rank determinations in the performed reductions are based on rank revealing QR-decompositions with column pivoting if fast = true
or the more reliable SVD-decompositions if fast = false
.
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 keyword argument atol3
specifies the absolute tolerance for the nonzero elements of B
and is only used if h2norm = false
for controllability tests of unstable eigenvalues. The keyword argument atolinf
is the absolute tolerance for the gain of G(λ)
at λ = ∞
. The used default value is atolinf = 0
. The default relative tolerance is n*ϵ
, 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
, atol2 = atol
and atol3 = atol
.
DescriptorSystems.gh2norm
— Functiongh2norm(sys, fast = true, offset = sqrt(ϵ), atol = 0, atol1 = atol, atol2 = atol, atolinf = atol, rtol = n*ϵ)
Compute for a descriptor system sys = (A-λE,B,C,D)
the H2
norm of its transfer function matrix G(λ)
. The H2
norm is infinite, if G(λ)
has unstable poles, or, for a continuous-time, the system has nonzero gain at infinity. If the pencil A-λE
has uncontrollable and/or unobservable unstable eigenvalues on the boundary of the stability domain, then a reduced order realization is determined first (see below) to eliminate these eigenvalues.
To check the stability, the eigenvalues of the pole pencil A-λE
must have real parts less than -β
for a continuous-time system or have moduli less than 1-β
for a discrete-time system, where β
is the stability domain boundary offset. The offset β
to be used can be specified via the keyword parameter offset = β
. The default value used for β
is sqrt(ϵ)
, where ϵ
is the working machine precision.
For a continuous-time system sys
with E
singular, a reduced order realization is determined first, without uncontrollable and unobservable finite and infinite eigenvalues of the pencil A-λE
. For a discrete-time system or for a system with invertible E
, a reduced order realization is determined first, without uncontrollable and unobservable finite eigenvalues of the pencil A-λE
. The rank determinations in the performed reductions are based on rank revealing QR-decompositions with column pivoting if fast = true
or the more reliable SVD-decompositions if fast = false
.
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 keyword argument atolinf
is the absolute tolerance for the gain of G(λ)
at λ = ∞
. The used default value is atolinf = 0
. The default relative tolerance is n*ϵ
, 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
.
DescriptorSystems.glinfnorm
— Functionglinfnorm(sys, hinfnorm = false, rtolinf = 0.001, fast = true, offset = sqrt(ϵ), atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ) -> (linfnorm, fpeak)
Compute for a descriptor system sys = (A-λE,B,C,D)
with the transfer function matrix G(λ)
the L∞
norm linfnorm
(i.e., the peak gain of G(λ)
) and the corresponding peak frequency fpeak
, where the peak gain is achieved. The L∞
norm is infinite if G(λ)
has poles on the stability domain boundary, i.e., on the extended imaginary axis, in the continuous-time case, or on the unit circle, in the discrete-time case. If the pencil A-λE
has uncontrollable and/or unobservable eigenvalues on the boundary of the stability domain, then a reduced order realization is determined first (see below) to eliminate these eigenvalues.
To check the lack of poles on the stability domain boundary, the eigenvalues of the pencil A-λE
must not have real parts in the interval [-β,β]
for a continuous-time system or must not have moduli within the interval [1-β,1+β]
for a discrete-time system, where β
is the stability domain boundary offset. The offset β
to be used can be specified via the keyword parameter offset = β
. The default value used for β
is sqrt(ϵ)
, where ϵ
is the working machine precision.
The keyword argument rtolinf
specifies the relative accuracy for the computed infinity norm. The default value used for rtolinf
is 0.001
.
If hinfnorm = true
, the H∞
norm is computed. In this case, the stability of the zeros of A-λE
is additionally checked and the H∞
norm is infinite for an unstable system. To check the stability, the eigenvalues of the pencil A-λE
must have real parts less than -β
for a continuous-time system or have moduli less than 1-β
for a discrete-time system.
For a continuous-time system sys
with E
singular, a reduced order realization is determined first, without uncontrollable and unobservable finite and infinite eigenvalues of the pencil A-λE
. For a discrete-time system or for a system with invertible E
, a reduced order realization is determined first, without uncontrollable and unobservable finite eigenvalues of the pencil A-λE
. The rank determinations in the performed reductions are based on rank revealing QR-decompositions with column pivoting if fast = true
or the more reliable SVD-decompositions if fast = false
.
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 n*ϵ
, 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
.
DescriptorSystems.ghinfnorm
— Functionghinfnorm(sys, rtolinf = 0.001, fast = true, offset = sqrt(ϵ), atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ) -> (hinfnorm, fpeak)
Compute for a descriptor system sys = (A-λE,B,C,D)
with the transfer function matrix G(λ)
the H∞
norm hinfnorm
(i.e., the peak gain of G(λ)
) and the corresponding peak frequency fpeak
, where the peak gain is achieved. The H∞
norm is infinite if G(λ)
has unstable poles. If the pencil A-λE
has uncontrollable and/or unobservable unstable eigenvalues, then a reduced order realization is determined first (see below) to eliminate these eigenvalues.
To check the stability, the eigenvalues of the pencil A-λE
must have real parts less than -β
for a continuous-time system or have moduli less than 1-β
for a discrete-time system, where β
is the stability domain boundary offset. The offset β
to be used can be specified via the keyword parameter offset = β
. The default value used for β
is sqrt(ϵ)
, where ϵ
is the working machine precision.
The keyword argument rtolinf
specifies the relative accuracy for the computed infinity norm. The default value used for rtolinf
is 0.001
.
For a continuous-time system sys
with E
singular, a reduced order realization is determined first, without uncontrollable and unobservable finite and infinite eigenvalues of the pencil A-λE
. For a discrete-time system or for a system with invertible E
, a reduced order realization is determined first, without uncontrollable and unobservable finite eigenvalues of the pencil A-λE
. The rank determinations in the performed reductions are based on rank revealing QR-decompositions with column pivoting if fast = true
or the more reliable SVD-decompositions if fast = false
.
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 n*ϵ
, 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
.
DescriptorSystems.gnugap
— Functiongnugap(sys1, sys2; freq = ω, rtolinf = 0.00001, fast = true, offset = sqrt(ϵ),
atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ) -> (nugapdist, fpeak)
Compute the ν-gap distance nugapdist
between two descriptor systems sys1 = (A1-λE1,B1,C1,D1)
and sys2 = (A2-λE2,B2,C2,D2)
and the corresponding frequency fpeak
(in rad/TimeUnit), where the ν-gap distance achieves its peak value.
If freq = missing
, the resulting nugapdist
satisfies 0 <= nugapdist <= 1
. The value nugapdist = 1
results, if the winding number is different of zero in which case fpeak = []
.
If freq = ω
, where ω
is a given vector of real frequency values, the resulting nugapdist
is a vector of pointwise ν-gap distances of the dimension of ω
, whose components satisfies 0 <= maximum(nugapdist) <= 1
. In this case, fpeak
is the frequency for which the pointwise distance achieves its peak value. All components of nugapdist
are set to 1 if the winding number is different of zero in which case fpeak = []
.
The stability boundary offset, β
, to be used to assess the finite zeros which belong to the boundary of the stability domain can be specified via the keyword parameter offset = β
. Accordingly, for a continuous-time system, these are the finite zeros having real parts within the interval [-β,β]
, while for a discrete-time system, these are the finite zeros having moduli within the interval [1-β,1+β]
. The default value used for β
is sqrt(ϵ)
, where ϵ
is the working machine precision.
Pencil reduction algorithms are employed to compute range and coimage spaces which perform rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true
or the more reliable SVD-decompositions if fast = false
.
The keyword arguments atol1
, atol2
and rtol
, specify, respectively, the absolute tolerance for the nonzero elements of A1
, A2
, B1
, B2
, C1
, C2
, D1
and D2
, the absolute tolerance for the nonzero elements of E1
and E2
, and the relative tolerance for the nonzero elements of all above matrices. The default relative tolerance is n*ϵ
, where ϵ
is the working machine epsilon and n
is the maximum of the orders of the systems sys1
and sys2
. The keyword argument atol
can be used to simultaneously set atol1 = atol
, atol2 = atol
.
The keyword argument rtolinf
specifies the relative accuracy to be used to compute the ν-gap as the infinity norm of the relevant system according to [1]. The default value used for rtolinf
is 0.00001
.
Method: The evaluation of ν-gap uses the definition proposed in [1], extended to generalized LTI (descriptor) systems. The computation of winding number is based on enhancements covering zeros on the boundary of the stability domain and infinite zeros.
References:
[1] G. Vinnicombe. Uncertainty and feedback: H∞ loop-shaping and the ν-gap metric. Imperial College Press, London, 2001.
DescriptorSystems.freqresp
— Function H = freqresp(sys, ω)
Compute the frequency response H
of the descriptor system sys = (A-λE,B,C,D)
at the real frequencies ω
.
For a system with p
outputs and m
inputs, and for N
frequency values in ω
, H
is a p × m × N
array, where H[:,:,i]
contains the i
-th value of the frequency response computed as C*((w[i]*E - A)^-1)*B + D
, where w[i] = im*ω[i]
for a continuous-time system and w[i] = exp(im*ω[i]*|Ts|)
for a discrete-time system with sampling time Ts
.
Method: For an efficient evaluation of C*((w[i]*E - A)^-1)*B + D
for many values of w[i]
, a preliminary orthogonal coordinate transformation is performed such that for the input-output equivalent transformed system sys = (At-λEt,Bt,Ct,Dt)
, the matrix w[i]*Et - At
is upper Hessenberg. This allows an efficient computation of the frequency response using the Hessenberg-form based approach proposed in [1].
Reference:
[1] Laub, A.J., "Efficient Multivariable Frequency Response Computations", IEEE Transactions on Automatic Control, AC-26 (1981), pp. 407-408.
DescriptorSystems.timeresp
— Functiontimeresp(sys, u, t, x0 = 0; interpolation = "zoh", state_history = false,
fast = true, atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ) -> (y, tout, x)
Compute the time response of a proper descriptor system sys = (A-λE,B,C,D)
to the input signals described by u
and t
. The time vector t
consists of regularly spaced time samples. The matrix u
has as many columns as the inputs of sys
and its i
-th row specifies the input values at time t[i]
. For discrete-time models, u
should be sampled at the same rate as sys
if sys.Ts > 0
and t
must have all time steps equal to sys.Ts
or can be set to an empty vector. The vector x0
specifies the initial state vector at time t[1]
and is set to zero when omitted.
The matrix y
contains the resulting time history of the outputs of sys
and the vector tout
contains the corresponding values of the time samples. The i
-th row of y
contains the output values at time tout[i]
. If the keyword parameter value state_history = false
is used, then the matrix x
contains the resulting time history of the state vector and its i
-th row contains the state values at time tout[i]
. By default, the state history is not computed and x = nothing
.
For continuous-time models, the input values are interpolated between samples. By default, zero-order hold based interpolation is used. The linear interpolation method can be selected using the keyword parameter interpolation = "foh"
.
By default, the uncontrollable infinite eigenvalues and simple infinite eigenvalues of the pair (A,E)
are eliminated. The underlying pencil manipulation algorithms employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true
(default), 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 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 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
.
DescriptorSystems.stepresp
— Functionstepresp(sys[, tfinal]; x0 = zeros(sys.nx), ustep = ones(sys.nu), timesteps = 100,
state_history = false, fast = true, atol = 0, atol1 = atol, atol2 = atol, rtol = n*ϵ) -> (y, tout, x)
Compute the step response of a proper descriptor system sys = (A-λE,B,C,D)
to step input signals. The final time tfinal
, if not specified, is set to 10 for a continuous-time system or to abs(sys.Ts)*timesteps
for a discrete-time system, where the keyword argument timesteps
specifies the number of desired simulation time steps (default: timesteps = 100
). The keyword argument ustep
is a vector with as many components as the inputs of sys
and specifies the desired amplitudes of step inputs (default: all components are set to 1). The keyword argument x0
is a vector which specifies the initial state vector at time 0
and is set to zero when omitted.
If ns
is the total number of simulation values, n
the number of state components, p
the number of system outputs and m
the number of system inputs, then the resulting ns×p×m
array y
contains the resulting time histories of the outputs of sys
, such that y[:,:,j]
is the time response for the j
-th input set to ustep[j]
and the rest of inputs set to zero. The vector tout
contains the corresponding values of the time samples. The i
-th row y[i,:,j]
contains the output values at time tout[i]
of the j
-th step response. If the keyword parameter value state_history = true
is used, then the resulting ns×n×m
arrayx
contains the resulting time histories of the state vector and the i
-th row x[i,:,j]
contains the state values at time tout[i]
of the j
-th step response. By default, the state history is not computed and x = nothing
.
By default, the uncontrollable infinite eigenvalues and simple infinite eigenvalues of the pair (A,E)
are eliminated. The underlying pencil manipulation algorithms employ rank determinations based on either the use of rank revealing QR-decomposition with column pivoting, if fast = true
(default), 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 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 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
.
DescriptorSystems.gbalqual
— Functionqs = gbalqual(sys; SysMat = false)
Compute the 1-norm based scaling quality of the matrices of the descriptor system sys = (A-λE,B,C)
.
If SysMat = false
, the resulting qs
is computed as
qs = max(qS(A),qS(E),qS(B),qS(C)) ,
where qS(⋅)
is the scaling quality measure defined in Definition 5.5 of [1] for nonnegative matrices, extended to also cover matrices with zero rows or columns.
If SysMat = true
, the resulting qs
is computed as
qs = qS(S) ,
where S
is the system matrix defined as
S = ( abs(A)+abs(E) abs(B) )
( abs(C) 0 )
A large value of qs
indicates a possibly poorly scaled state-space model. For a standard system with E = I
, the above formulas are used assuming E = 0
.
[1] F.M.Dopico, M.C.Quintana and P. van Dooren, "Diagonal scalings for the eigenstructure of arbitrary pencils", SIMAX, 43:1213-1237, 2022.