Manipulation of structured linear matrix pencils

FunctionDescription
sreduceBFReduction to the basic condensed form [B A-λE; D C] with E upper triangular and nonsingular.
sklfComputation of the Kronecker-like form exhibiting the full Kronecker structure
sklf_rightComputation of the Kronecker-like form exhibiting the right Kronecker structure
sklf_leftComputation of the Kronecker-like form exhibiting the left Kronecker structure
gsklfComputation of several row partition preserving special Kronecker-like forms
MatrixPencils.sklfFunction
sklf(A, E, B, C, D; fast = true, finite_infinite = false, ut = false, atol1 = 0, atol2 = 0, rtol, withQ = true, withZ = true) -> (F, G, Q, Z, νr, μr, νi, nf, νl, μl)

Reduce the structured matrix pencil M - λN

           | A-λE | B | 
  M - λN = |------|---|
           |  C   | D |

to an equivalent form F - λG = Q'*(M - λN)*Z using orthogonal or unitary transformation matrices Q and Z such that the transformed matrices F and G are in the following Kronecker-like form exhibiting the complete Kronecker structure:

            |  Mr-λNr  |     *      |    *    |
            |----------|------------|---------|
  F - λG =  |    O     | Mreg-λNreg |    *    |
            |----------|------------|---------|
            |    O     |     0      |  Ml-λNl |

The full row rank pencil Mr-λNr is in a staircase form, contains the right Kronecker indices of the pencil M-λNand has the form

     Mr-λNr  = | Br | Ar-λEr |,

where Er is upper triangular and nonsingular. The nr-dimensional vectors νr and μr contain the row and, respectively, column dimensions of the blocks of the staircase form Mr-λNr such that the i-th block has dimensions νr[i] x μr[i] and has full row rank. The difference μr[i]-νr[i] for i = 1, 2, ..., nr is the number of elementary Kronecker blocks of size (i-1) x i. If ut = true, the full row rank diagonal blocks of Mr are reduced to the form [0 X] with X upper triangular and nonsingular and the full column rank supradiagonal blocks of Nr are reduced to the form [Y; 0] with Y upper triangular and nonsingular.

The pencil Mreg-λNreg is regular, contains the infinite and finite elementary divisors of M-λN and has the form

                   | Mi-λNi |    *    |
     Mreg-λNreg  = |--------|---------|, if `finite_infinite = false`, or the form
                   |   0    |  Mf-λNf |


                   | Mf-λNf |    *    |
     Mreg-λNreg  = |--------|---------|, if `finite_infinite = true`, 
                   |   0    |  Mi-λNi |

where: (1) Mi-λNi, in staircase form, contains the infinite elementary divisors of M-λN, Mi is upper triangular if ut = true and nonsingular, and Ni is upper triangular and nilpotent; (2) Mf-λNf contains the infinite elementary divisors of M-λN and Nf is upper triangular and nonsingular. The ni-dimensional vector νi contains the dimensions of the square blocks of the staircase form Mi-λNi such that the i-th block has dimensions νi[i] x νi[i]. If finite_infinite = true, the difference νi[i]-νi[i+1] for i = 1, 2, ..., ni is the number of infinite elementary divisors of degree i (with νi[ni] = 0). If finite_infinite = false, the difference νi[ni-i+1]-νi[ni-i] for i = 1, 2, ..., ni is the number of infinite elementary divisors of degree i (with νi[0] = 0).

The full column rank pencil Ml-λNl, in a staircase form, contains the left Kronecker indices of M-λN and has the form

               | Al-λEl |
     Ml-λNl  = |--------|,
               |   Cl   |

where El is upper triangular and nonsingular. The nl-dimensional vectors νl and μl contain the row and, respectively, column dimensions of the blocks of the staircase form Ml-λNl such that the j-th block has dimensions νl[j] x μl[j] and has full column rank. The difference νl[nl-j+1]-μl[nl-j+1] for j = 1, 2, ..., nl is the number of elementary Kronecker blocks of size j x (j-1). If ut = true, the full column rank diagonal blocks of Ml are reduced to the form [X; 0] with X upper triangular and nonsingular and the full row rank supradiagonal blocks of Nl are reduced to the form [0 Y] with Y upper triangular and nonsingular.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of M, the absolute tolerance for the nonzero elements of N, and the relative tolerance for the nonzero elements of M and N. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed left orthogonal or unitary transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed right orthogonal or unitary transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

source
MatrixPencils.sklf_rightFunction
sklf_right(A, E, B, C, D; fast = true, ut = false, atol1 = 0, atol2 = 0, rtol, withQ = true, withZ = true) -> (F, G, Q, Z, νr, μr, nf, ν, μ)

Reduce the structured matrix pencil M - λN

           | A-λE | B | 
  M - λN = |------|---|
           |  C   | D |

to an equivalent form F - λG = Q'*(M - λN)*Z using orthogonal or unitary transformation matrices Q and Z such that the transformed matrices F and G are in the following Kronecker-like form exhibiting the right and finite Kronecker structures:

            |  Mr-λNr    |     *      |    *     |
            |------------|------------|----------|
  F - λG =  |    O       |   Mf-λNf   |    *     |
            |------------|------------|----------|
            |    O       |     0      | Mil-λNil |

The full row rank pencil Mr-λNr, in a staircase form, contains the right Kronecker indices of M-λN and has the form

     Mr-λNr  = | Br | Ar-λEr |,

where Er is upper triangular and nonsingular. The nr-dimensional vectors νr and μr contain the row and, respectively, column dimensions of the blocks of the staircase form Mr-λNr such that i-th block has dimensions νr[i] x μr[i] and has full row rank. The difference μr[i]-νr[i] for i = 1, 2, ..., nr is the number of elementary Kronecker blocks of size (i-1) x i. If ut = true, the full row rank diagonal blocks of Mr are reduced to the form [0 X] with X upper triangular and nonsingular and the full column rank supradiagonal blocks of Nr are reduced to the form [Y; 0] with Y upper triangular and nonsingular.

The nf x nf pencil Mf-λNf is regular and contains the finite elementary divisors of M-λN. Nf is upper triangular and nonsingular.

The full column rank pencil Mil-λNil is in a staircase form, and contains the left Kronecker indices and infinite elementary divisors of the pencil M-λN. The nb-dimensional vectors ν and μ contain the row and, respectively, column dimensions of the blocks of the staircase form Mil-λNil such that the i-th block has dimensions ν[i] x μ[i] and has full column rank. The difference ν[nb-j+1]-μ[nb-j+1] for j = 1, 2, ..., nb is the number of elementary Kronecker blocks of size j x (j-1). The difference μ[nb-j+1]-ν[nb-j] for j = 1, 2, ..., nb is the number of infinite elementary divisors of degree j (with ν[0] = 0). If ut = true, the full column rank diagonal blocks of Mil are reduced to the form [X; 0] with X upper triangular and nonsingular and the full row rank supradiagonal blocks of Nil are reduced to the form [0 Y] with Y upper triangular and nonsingular.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of M, the absolute tolerance for the nonzero elements of N, and the relative tolerance for the nonzero elements of M and N. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed left orthogonal or unitary transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed right orthogonal or unitary or unitary transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

Note: If the pencil M - λN has full row rank, then the regular pencil Mil-λNil is in a staircase form with square upper triangular diagonal blocks (i.e.,μ[j] = ν[j]), and the difference ν[nb-j+1]-ν[nb-j] for j = 1, 2, ..., nb is the number of infinite elementary divisors of degree j (with ν[0] = 0).

source
MatrixPencils.sklf_leftFunction
sklf_left(A, E, B, C, D; fast = true, ut = false, atol1 = 0, atol2 = 0, rtol, withQ = true, withZ = true) -> (F, G, Q, Z, ν, μ, nf, νl, μl)

Reduce the structured matrix pencil M - λN

           | A-λE | B | 
  M - λN = |------|---|
           |  C   | D |

to an equivalent form F - λG = Q'*(M - λN)*Z using orthogonal or unitary transformation matrices Q and Z such that the transformed matrices F and G are in the following Kronecker-like form exhibiting the finite and left Kronecker structures:

            |  Mri-λNri  |     *      |    *    |
            |------------|------------|---------|
  F - λG =  |    O       |   Mf-λNf   |    *    |
            |------------|------------|---------|
            |    O       |     0      |  Ml-λNl |

The full row rank pencil Mri-λNri is in a staircase form, and contains the right Kronecker indices and infinite elementary divisors of the pencil M-λN.

The nf x nf pencil Mf-λNf is regular and contains the finite elementary divisors of M-λN. Nf is upper triangular and nonsingular.

The full column rank pencil Ml-λNl, in a staircase form, contains the left Kronecker indices of M-λN and has the form

               | Al-λEl |
     Ml-λNl  = |--------|,
               |   Cl   |

where El is upper triangular and nonsingular.

The nb-dimensional vectors ν and μ contain the row and, respectively, column dimensions of the blocks of the staircase form Mri-λNri such that i-th block has dimensions ν[i] x μ[i] and has full row rank. The difference μ[i]-ν[i] for i = 1, 2, ..., nb is the number of elementary Kronecker blocks of size (i-1) x i. The difference ν[i]-μ[i+1] for i = 1, 2, ..., nb is the number of infinite elementary divisors of degree i (with μ[nb+1] = 0). If ut = true, the full row rank diagonal blocks of Mri are reduced to the form [0 X] with X upper triangular and nonsingular and the full column rank supradiagonal blocks of Nri are reduced to the form [Y; 0] with Y upper triangular and nonsingular.

The nl-dimensional vectors νl and μl contain the row and, respectively, column dimensions of the blocks of the staircase form Ml-λNl such that the j-th block has dimensions νl[j] x μl[j] and has full column rank. The difference νl[nl-j+1]-μl[nl-j+1] for j = 1, 2, ..., nl is the number of elementary Kronecker blocks of size j x (j-1).

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of M, the absolute tolerance for the nonzero elements of N, and the relative tolerance for the nonzero elements of M and N. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed left orthogonal or unitary transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed right orthogonal or unitary transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

Note: If the pencil M - λN has full column rank, then the regular pencil Mri-λNri is in a staircase form with square upper triangular diagonal blocks (i.e.,μ[i] = ν[i]), and the difference ν[i+1]-ν[i] for i = 1, 2, ..., nb is the number of infinite elementary divisors of degree i (with ν[nb+1] = 0).

source
MatrixPencils.sreduceBFFunction
sreduceBF(A, E, B, C, D; fast = true, atol = 0, rtol, withQ = true, withZ = true) -> F, G, Q, Z, n, m, p

Reduce the partitioned matrix pencil M - λN

           | A-λE | B | ndx
  M - λN = |------|---|
           |  C   | D | ny
              nx    nu

to an equivalent basic form F - λG = Q'*(M - λN)*Z using orthogonal transformation matrices Q and Z such that M - λN is transformed into the following standard form

           | B1 | A1-λE1 | n
  F - λG = |----|--------|   ,        
           | D1 |   C1   | p
             m      n

where E1 is an nxn non-singular matrix, and A1, B1, C1, D1 are nxn-, nxm-, pxn- and pxm-dimensional matrices, respectively. The order n of E1 is equal to the numerical rank of E determined using the absolute tolerance atol and relative tolerance rtol. The dimensions m and p are computed as m = nu + (nx-n) and p = ny + (ndx-n).

If fast = true, E1 is determined upper triangular using a rank revealing QR-decomposition with column pivoting of E and n is evaluated as the number of nonzero diagonal elements of the R factor, whose magnitudes are greater than tol = max(atol,abs(R[1,1])*rtol). If fast = false, E1 is determined diagonal using a rank revealing SVD-decomposition of E and n is evaluated as the number of singular values greater than tol = max(atol,smax*rtol), where smax is the largest singular value. The rank decision based on the SVD-decomposition is generally more reliable, but the involved computational effort is higher.

The performed left orthogonal or unitary transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed right orthogonal or unitary transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

source
MatrixPencils.gsklfFunction
gsklf(A, E, B, C, D; disc = false, jobopt = "none", offset = sqrt(ϵ), fast = true, atol1 = 0, atol2 = 0, rtol, 
     withQ = true, withZ = true) -> (F, G, Q, Z, dimsc, nmszer, nizer)

Reduce the structured matrix pencil M - λN

           | A-λE | B | 
  M - λN = |------|---|
           |  C   | D |

with the n x n pencil A-λE regular and [E B] of full row rank n, to an equivalent form F - λG = diag(Q',I)*(M - λN)*Z using orthogonal or unitary transformation matrices Q and Z such that the transformed matrices F and G are in the following special, row partition preserving, Kronecker-like form

            | Ar-λEr    *     *   *   |
            |   0     Al-λEl  Bl  *   |
   F - λG = |   0       0     0   Bn  | ,    
            |-------------------------|
            |   0       Cl    Dl  *   |

where: (1) Ar-λEr is a mr x nr full row rank pencil, which contains the right Kronecker structure of M - λN and that part of its zeros which lie outside of the domain of interest Cb of the complex plane defined by the keyword arguments disc and jobopt; (2) Al-lambda*El is an nl x nl regular square pencil with El invertible and upper triangular; (3) Bl is an nl x ml matrix, where ml is the normal rank of the transfer function matrix C*inv(λE-A)*B+D; (4) Bn is an invertible matrix of order nsinf.

The column dimensions (nr, nl, ml, nsinf) are returned in dimsc.

The full column rank subpencil

            | Al-λEl | Bl | 
Ml - λNl := |--------|----|
            |   Cl   | Dl |

contains the left Kronecker structure of M - λN and its zeros lying in the domain of interest Cb, and has the property that rank [Al-λEl Bl] = nl provided rank [A-λE B] = n for all finite λ ∈ Cb (finite stabilizability with respect to Cb). The number of finite zeros lying on the boundary of Cb, if available, is returned in nmszer, otherwise nmszer is set to missing. If disc = false, the number of infinite zeros, if available, is returned in nizer, otherwise nizer is set to missing. If disc = true, nizer is set to 0. The boundary offset β to be used to assess the zeros on the boundary of Cb is specified by the keyword parameter offset. Accordingly, if disc = false, then the boundary of Cb contains the complex numbers with real parts within the interval [-β,β], while if disc = true, then the boundary of Cb contains the complex numbers with moduli within the interval [1-β,1+β]. The default value used for β is sqrt(ϵ) (the machine precision of the element type of A).

The domain of interest Cb for the zeros of the pencil Ml - λNl is defined as follows:

if jobopt = none, then Cb is the empty set and no zeros of M - λN are included in Ml - λNl;

if jobopt = all, then Cb is the extended complex plane (including infinity) and all zeros of M - λN are included in Ml - λNl;

if jobopt = finite, then Cb is the complex plane (without infinity) and all finite zeros of M - λN are included in Ml - λNl;

if jobopt = infinite, then Cb is the point at infinity and all infinite zeros of M - λN are included in Ml - λNl;

if jobopt = stable, then, for disc = false, Cb is the set of complex numbers with real parts less than , while for disc = true, Cb is the set of complex numbers with moduli less than 1-β and all finite zeros of M - λN in Cb are included in Ml - λNl;

if jobopt = unstable, then, for disc = false, Cb is the set of complex numbers with real parts greater than β or infinite, while for disc = true, Cb is the set of complex numbers with moduli greater than 1+β or infinite and all finite and infinite zeros of M - λN in Cb are included in Ml - λNl;

if jobopt = s-unstable, then, for disc = false, Cb is the set of complex numbers with real parts greater than or infinite, while for disc = true, Cb is the set of complex numbers with moduli greater than 1-β or infinite and all finite and infinite zeros of M - λN in Cb are included in Ml - λNl.

The keyword arguments atol1, atol2, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of M, the absolute tolerance for the nonzero elements of N, and the relative tolerance for the nonzero elements of M and N. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed left orthogonal or unitary transformations are accumulated in the n x n matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed right orthogonal or unitary transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

Method: The reduction algorithm of [1] has been adapted to deal with several zero selection options.

References

[1] Oara, C. Constructive solutions to spectral and inner–outer factorizations with respect to the disk. Automatica, 41:1855–1866, 2005.

source
MatrixPencils.sklf_right!Function
sklf_right!(A, E, B, F, C, G, D, H; fast = true, atol1 = 0, atol2 = 0, atol3 = 0, rtol, withQ = true, withZ = true) -> (Q, Z, nc)

Reduce the partitioned full row rank matrix pencil

  M - λN = | B-λF | A-λE | n
              m      n

with A-λE regular, to an equivalent form Mt - λNt = Q'*(M - λN)*Z = | Bt-λFt | At-λEt | using orthogonal or unitary transformation matrices Q and Z, such that M - λN is transformed into a Kronecker-like form | Bt-λFt | At-λEt | exhibiting its right, finite and infinite Kronecker structures (also known as the strong controllability form):

                         | Bc-λFc | Ac-λEc     *     | nc
  | Bt-λFt | At-λEt | =  |--------|------------------|
                         |   0    |   0     Auc-λEuc | n-nc
                             m        nc      n-nc

The matrices Bt, Ft, At, Et, determined from [Bt At] = Q'*[B A]*Z, [Ft Et] = Q'*[F E]*Z, are returned in B, F, A and E, respectively. Furthermore, the matrices Ct, Dt, Gt, Ht, determined from the compatibly partitioned matrices [Dt Ct] = [D C]*Z and [Ht Gt] = [D C]*Z, are returned in C, D, G and H, respectively (unless C, D, G and H are set to missing).

The subpencil | Bc-λFc | Ac-λEc | has full row rank nc and the (n-nc) x (n-nc) subpencil Auc-λEuc contains the finite and infinite eigenvalues of M - λN (also called the uncontrollable eigenvalues of M - λN).

The (m+n) x (m+n) orthogonal matrix Z has the partitioned form

       | Z11 |  0  |  *  | m 
   Z = |-----|-----|-----|
       |  *  |  *  |  *  | n
          m    nc    n-nc

with the leading m x m submatrix Z11 invertible and upper triangular.

Note: If Ct-λGt is partitioned as [ Cc-λGc | Cuc-λGuc ], then the following relation is fulfilled

  (Cc-λGc)*inv(λEc-Ac)(Bc-λFc) + Dt-λHt = ((C-λG)*inv(λE-A)(B-λF) + D-λH)*Z11

and the structured pencil linearization (Ac-λEc,Bc-λFc,Cc-λGc,Dt-λHt) is strongly controllable.

The keyword arguments atol1, atol2, , atol3, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of A and B, the absolute tolerance for the nonzero elements of E and F, and the relative tolerance for the nonzero elements of A, B, E and F. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary left transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed orthogonal or unitary right transformations are accumulated in the matrix Z if withZ = true. If withZ = false, Z contains the upper triangular matrix Z11.

source
sklf_right!(A, E, B, C; fast = true, atol1 = 0, atol2 = 0, atol3 = 0, rtol, withQ = true, withZ = true) -> (Q, Z, νr, nc, nfuc, niuc)

Reduce the partitioned full row rank matrix pencil

  M - λN = | B | A-λE | n
             m    n

with A-λE regular, to an equivalent form F - λG = Q'*(M - λN)*diag(I,Z) using orthogonal or unitary transformation matrices Q and Z such that M - λN is transformed into a Kronecker-like form | Bt | At-λEt | exhibiting its right, finite and infinite Kronecker structures (also known as the generalized controllability staircase form):

                     |  Bc | Ac-λEc     *          *      | nc
  | Bt | At-λEt | =  |-----|------------------------------|
                     |  0  |  0     Afuc-λEfuc     *      | nfuc
                     |  0  |  0         0      Aiuc-λEiuc | niuc
                        m     nc       nfuc        niuc

Bt = Q'*B, At = Q'*A*Zand Et = Q'*E*Z are returned in B, A and E, respectively, and C*Z is returned in C (unless C = missing).

The subpencil | Bc | Ac-λEc | has full row rank nc, is in a staircase form, and contains the right Kronecker indices of M - λN. The nr-dimensional vector νr contains the row and column dimensions of the blocks of the staircase form | Bc | Ac-λI | such that i-th block has dimensions νr[i] x νr[i-1] (with νr[0] = m) and has full row rank. The difference νr[i-1]-νr[i] for i = 1, 2, ..., nr is the number of elementary Kronecker blocks of size (i-1) x i.

The nfuc x nfuc subpencil Afuc-λEfuc contains the finite eigenvalues of M - λN (also called the uncontrollable finite eigenvalues of A - λE).

The niuc x niuc subpencil Aiuc-λEiuc contains the infinite eigenvalues of M - λN (also called the uncontrollable infinite eigenvalues of A - λE).

The keyword arguments atol1, atol2, atol3, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of A, the absolute tolerance for the nonzero elements of E, the absolute tolerance for the nonzero elements of B, and the relative tolerance for the nonzero elements of A, E and B. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary left transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed orthogonal or unitary right transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

source
sklf_right!(A, B, C; fast = true, atol1 = 0, atol2 = 0, rtol, withQ = true) -> (Q, νr, nc, nuc)

Reduce the partitioned full row rank matrix pencil

  M - λN = | B | A-λI | n
             m    n

to an equivalent form F - λG = Q'*(M - λN)*diag(I,Q) using orthogonal or unitary transformation matrix Q such that M - λN is transformed into a Kronecker-like form | Bt | At-λI | exhibiting its right and finite Kronecker structures (also known as the controllability staircase form):

                    |  Bc | Ac-λI     *    | nc
  | Bt | At-λI | =  |-----|----------------|
                    |  0  |  0     Auc-λI  | nuc
                       m     nc      nuc

Bt = Q'*B and At = Q'*A*Q are returned in B and A, respectively, and C*Q is returned in C (unless C = missing).

The subpencil | Bc | Ac-λI | has full row rank nc, is in a staircase form, and contains the right Kronecker indices of M - λN. The nr-dimensional vector νr contains the row and column dimensions of the blocks of the staircase form | Bc | Ac-λI | such that i-th block has dimensions νr[i] x νr[i-1] (with νr[0] = m) and has full row rank. The difference νr[i-1]-νr[i] for i = 1, 2, ..., nr is the number of elementary Kronecker blocks of size (i-1) x i.

The nuc x nuc matrix Auc contains the (finite) eigenvalues of M - λN (also called the uncontrollable eigenvalues of A).

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 B, and the relative tolerance for the nonzero elements of A and B. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing.

source
MatrixPencils.sklf_right2!Function
sklf_right2!(A, B, m1, C; fast = true, atol1 = 0, atol2 = 0, rtol, withQ = true) -> (Q, νr, nc, nuc)

Reduce the partitioned full row rank matrix pencil

  M - λN = [ B | A-λI ] n  = [ B1  B2 | A-λI ] n
             m    n            m1 m-m1   n

to an equivalent form F - λG = Q'*(M - λN)*diag(I,Q) using orthogonal or unitary transformation matrix Q such that M - λN is transformed into a Kronecker-like form [ Bt | At-λI ] exhibiting its right and finite Kronecker structures (also known as the controllability staircase form):

                    [  Bc | Ac-λI     *    ] nc
  [ Bt | At-λI ] =  [-----|----------------]
                    [  0  |  0     Auc-λI  ] nuc
                       m     nc      nuc

Bt = Q'*B and At = Q'*A*Q are returned in B and A, respectively, and C*Q is returned in C (unless C = missing).

The subpencil [ Bc | Ac-λI ] has full row rank nc with [ Bc | Ac ] = [ Bc1 Bc2 | Ac] in a staircase form

                    m1  m-m1  νr[1] νr[2] . . .  νr[p-2] 
                  [ B11 B12 | A11   A12   . . .  A1,p-2  A1,p-1  A1p ]  νr[1]
                  [  0  B22 | A21   A22   . . .  A2,p-2  A2,p-1  A2p ]  νr[2]
                  [  0   0  | A31   A32   . . .  A3,p-2  A3,p-1  A3p ]  νr[3]
                  [  0   0  |  0    A42   . . .  A4,p-2  A4,p-1  A4p ]  νr[4]
[ Bc1 Bc2 | Ac] = [  0   0  |  .     .    . . .    .       .      .  ]   .
                  [  0   0  |  .     .      . .    .       .      .  ]   .
                  [  0   0  |  .     .        .    .       .      .  ]   .
                  [  0   0  |  0     0    . . .  Ap,p-2  Ap,p-1  App ]  νr[p]

where the blocks B11, B22, A31, ..., Ap,p-2 have full row ranks. The p-dimensional vector νr, with p = 2nr even, contains the row dimensions of the blocks. The difference νr[2i-1]+νr[2i]-νr[2i+1]-νr[2i+2] for i = 1, 2, ..., nr is the number of elementary Kronecker blocks of size (i-1) x i.

The nuc x nuc matrix Auc contains the (finite) eigenvalues of M - λN (also called the uncontrollable eigenvalues of A).

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 B, and the relative tolerance for the nonzero elements of A and B. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing.

Method: The implemented algorithm [1] represents a specialization of the controllability staircase algorithm of [2] to the special structure of the input matrix B = [B1,B2].

References:

[1] Varga, A. Reliable algorithms for computing minimal dynamic covers. Proc. CDC'2003, Hawaii, 2003.

[2] Varga, A. Numerically stable algorithm for standard controllability form determination. Electronics Letters, vol. 17, pp. 74-75, 1981.

source
MatrixPencils.sklf_left!Function
sklf_left!(A, E, C, G, B, F, D, H; fast = true, atol1 = 0, atol2 = 0, atol3 = 0, rtol, withQ = true, withZ = true) -> (Q, Z, no)

Reduce the partitioned full row rank matrix pencil

M - λN = | A-λE | n
         | C-λG | p
            n

with A-λE regular, to an equivalent form F - λG = Q'*(M - λN)*Z using orthogonal or unitary transformation matrices Q and Z such that M - λN is transformed into a staircase form exhibiting the left, finite and infinite Kronecker structures (also known as the strong observability form):

                | Auo-λEuo  |   *    | n-no
  | At-λEt |    |   0       | Ao-λEo | no
  |--------| =  |-----------|--------|
  | Ct-λGt |    |   0       | Co-λGo | p
                   n-no         no

The matrices Ct, Gt, At, Et, determined from [At; Ct] = Q'*[A; C]*Z, [Et; Gt] = Q'*[E; G]*Z, are returned in C, G, A and E, respectively. Furthermore, the matrices Bt, Dt, Ft, Ht determined from the compatibly partitioned [Bt; Dt] = Q'*[B; D] and [Ft; Ht] = Q'*[F; H], are returned in B, D, F, and H, respectively (unless B, D, F and H are set to missing).

The subpencil

   | Ao-λEo |   
   |--------|
   | Co-λGo |

has full column rank no and the (n-no) x (n-no) subpencil Auo-λEuo contains the finite and infinite eigenvalues of M - λN (also called the unobservable eigenvalues of M - λN).

The (n+p) x (n+p) orthogonal matrix Q has the partitioned form

       |  *  |  *  |  *   | n 
   Q = |-----|-----|------|
       |  *  |  0  | Q22  | p
         n-no  no     p

with the p x p trailing submatrix Q22 invertible and upper triangular.

Note: If Bt-λFt is partitioned as [ Buo-λFuo | Bo-λFo ], then the following relation is fulfilled

  (Co-λGo)*inv(λEo-Ao)(Bo-λFo) + Dt-λHt = Q22'*((C-λG)*inv(λE-A)(B-λF) + D-λH)

and the structured pencil linearization (Ao-λEo,Bo-λFo,Co-λGo,Dt-λHt) is strongly observable.

The keyword arguments atol1, atol2, , atol3, and rtol, specify, respectively, the absolute tolerance for the nonzero elements of A and C, the absolute tolerance for the nonzero elements of E and G, and the relative tolerance for the nonzero elements of A, C, E and G. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary left transformations are accumulated in the matrix Q if withQ = true. If withQ = false, Q contains the upper triangular matrix Q22. The performed orthogonal or unitary right transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

source
sklf_left!(A, E, C, B; fast = true, atol1 = 0, atol2 = 0, atol3 = 0, rtol, withQ = true, withZ = true) -> (Q, Z, μl, no, nfuo, niuo)

Reduce the partitioned full column rank matrix pencil

  M - λN = | A-λE | n
           |  C   | p
              n

to an equivalent form F - λG = diag(Q',I)*(M - λN)*Z using orthogonal or unitary transformation matrices Q and Z such that M - λN is transformed into a Kronecker-like form exhibiting its infinite, finite and left Kronecker structures (also known as the generalized observability staircase form):

                | Aiuo-λEiuo     *       |   *    | niuo
                |    0       Afuo-λEfuo  |   *    | nfuo
  | At-λEt |    |    0           0       | Ao-λEo | no
  |--------| =  |------------------------|--------|
  |  Ct    |    |    0           0       |   Co   | p
                    niuo        nfuo         no

Ct = C*Z, At = Q'*A*Z and Et = Q'*E*Z are returned in C, A and E, respectively, and Q'*B is returned in B (unless B = missing).

The subpencil | Ao-λEo | has full column rank no, is in a staircase form, and contains the left Kronecker indices of M - λN. | Co | The nl-dimensional vector μl contains the row and column dimensions of the blocks of the staircase form such that i-th block has dimensions μl[nl-i] x μl[nl-i+1] (with μl[0] = p) and has full column rank. The difference μl[nl-i]-μl[nl-i+1] for i = 1, 2, ..., nl is the number of elementary Kronecker blocks of size i x (i-1).

The niuo x niuo subpencil Aiuo-λEiuo contains the infinite eigenvalues of M - λN (also called the unobservable infinite eigenvalues of A-λE).

The nfuo x nfuo subpencil Afuo-λEfuo contains the finite eigenvalues of M - λN (also called the unobservable finite eigenvalues of A-λE).

The keyword arguments atol1, atol2, atol3 and rtol, specify, respectively, the absolute tolerance for the nonzero elements of A, the absolute tolerance for the nonzero elements of E, the absolute tolerance for the nonzero elements of C, and the relative tolerance for the nonzero elements of A, E and C. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary left transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed orthogonal or unitary right transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

source
sklf_left!(A, C, B; fast = true, atol1 = 0, atol2 = 0, rtol, withQ = true) -> (Q, μl, no, nuo)

Reduce the partitioned full column rank matrix pencil

  M - λN = | A-λI | n
           |  C   | p
              n

to an equivalent form F - λL = diag(Q',I)*(M - λN)*Q using orthogonal or unitary transformation matrix Q such that M - λN is transformed into a Kronecker-like form exhibiting its finite and left Kronecker structures (also known as the observability staircase form):

               | Auo-λI  |  *    | nuo
  | At-λI |    |    0    | Ao-λI | no
  |-------| =  |---------|-------|
  |  Ct   |    |    0    |   Co  | p
                   nuo       no

Ct = C*Q and At = Q'*A*Q are returned in C and A, respectively, and Q'*B is returned in B (unless B = missing).

The subpencil | Ao-λI; Co | has full column rank no, is in a staircase form, and contains the left Kronecker indices of M - λN.

The nl-dimensional vector μl contains the row and column dimensions of the blocks of the staircase form such that i-th block has dimensions μl[nl-i] x μl[nl-i+1] (with μl[0] = p) and has full column rank. The difference μl[nl-i]-μl[nl-i+1] for i = 1, 2, ..., nl is the number of elementary Kronecker blocks of size i x (i-1).

The nuo x nuo matrix Auo contains the (finite) eigenvalues of M - λN (also called the unobservable eigenvalues of A).

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 C, and the relative tolerance for the nonzero elements of A and C. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing.

source
MatrixPencils.sklf_rightfin!Function
sklf_rightfin!(A, E, B, C; fast = true, atol1 = 0, atol2 = 0,  
               rtol, withQ = true, withZ = true) -> (Q, Z, νr, nc, nuc)

Reduce the partitioned full row rank matrix pencil

  M - λN = | B | A-λE | n
             m    n

with A-λE regular, to an equivalent form F - λG = Q'*(M - λN)*diag(I,Z) using orthogonal or unitary transformation matrices Q and Z such that M - λN is transformed into a staircase form | Bt | At-λEt | exhibiting the separation of its finite eigenvalues:

                     |  Bc | Ac-λEc     *    | nc
  | Bt | At-λEt | =  |-----|-----------------|
                     |  0  |  0     Auc-λEuc | nuc
                        m     nc      nuc

Bt = Q'*B, At = Q'*A*Zand Et = Q'*E*Z are returned in B, A and E, respectively, and C*Z is returned in C (unless C = missing). The resulting Et is upper triangular. If E is already upper triangular, then the preliminary reduction of E to upper triangular form is not performed.

The subpencil | Bc | Ac-λEc | has full row rank nc for all finite values of λ, is in a staircase form, and, if E is invertible, contains the right Kronecker indices of M - λN. The nr-dimensional vector νr contains the row and column dimensions of the blocks of the staircase form | Bc | Ac-λI | such that i-th block has dimensions νr[i] x νr[i-1] (with νr[0] = m) and has full row rank. If E is invertible, the difference νr[i-1]-νr[i] for i = 1, 2, ..., nr is the number of elementary Kronecker blocks of size (i-1) x i.

The nuc x nuc subpencil Auc-λEuc contains the finite eigenvalues of M - λN (also called the uncontrollable finite eigenvalues of A - λE). If E is singular, Auc-λEuc may also contain a part of the infinite eigenvalues of M - λN (also called the uncontrollable infinite eigenvalues of A - λE).

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 B, and the relative tolerance for the nonzero elements of A and B. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary left transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed orthogonal or unitary right transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

Note: This function, called with reversed input parameters E and A (i.e., instead A and E), performs the separation all infinite and nonzero finite eigenvalues of the pencil M - λN.

source
MatrixPencils.sklf_rightfin2!Function
sklf_rightfin2!(A, E, B, m1, C; fast = true, atol1 = 0, atol2 = 0,  
               rtol, withQ = true, withZ = true) -> (Q, Z, νr, nc, nuc)

Reduce the partitioned full row rank matrix pencil

  M - λN = [ B | A-λE ] n  = [ B1  B2 | A-λE ] n
             m    n            m1 m-m1   n

with A-λE regular, to an equivalent form F - λG = Q'*(M - λN)*diag(I,Z) using orthogonal or unitary transformation matrices Q and Z such that M - λN is transformed into a staircase form [ Bt | At-λEt ] exhibiting the separation of its finite eigenvalues:

                     [  Bc | Ac-λEc     *    ] nc
  [ Bt | At-λEt ] =  [-----|-----------------]
                     [  0  |  0     Auc-λEuc ] nuc
                        m     nc      nuc

Bt = Q'*B, At = Q'*A*Zand Et = Q'*E*Z are returned in B, A and E, respectively, and C*Z is returned in C (unless C = missing). The resulting Et is upper triangular. If E is already upper triangular, then the preliminary reduction of E to upper triangular form is not performed.

The subpencil [ Bc | Ac-λEc ] has full row rank nc for all finite values of λ, with [ Bc | Ac ] = [ Bc1 Bc2 | Ac] in a staircase form

                    m1  m-m1  νr[1] νr[2] . . .  νr[p-2] 
                  [ B11 B12 | A11   A12   . . .  A1,p-2  A1,p-1  A1p ]  νr[1]
                  [  0  B22 | A21   A22   . . .  A2,p-2  A2,p-1  A2p ]  νr[2]
                  [  0   0  | A31   A32   . . .  A3,p-2  A3,p-1  A3p ]  νr[3]
                  [  0   0  |  0    A42   . . .  A4,p-2  A4,p-1  A4p ]  νr[4]
[ Bc1 Bc2 | Ac] = [  0   0  |  .     .    . . .    .       .      .  ]   .
                  [  0   0  |  .     .      . .    .       .      .  ]   .
                  [  0   0  |  .     .        .    .       .      .  ]   .
                  [  0   0  |  0     0    . . .  Ap,p-2  Ap,p-1  App ]  νr[p]

where the blocks B11, B22, A31, ..., Ap,p-2 have full row ranks. The p-dimensional vector νr, with p = 2nr even, contains the row dimensions of the blocks. The difference νr[2i-1]+νr[2i]-νr[2i+1]-νr[2i+2] for i = 1, 2, ..., nr is the number of elementary Kronecker blocks of size (i-1) x i.

The nuc x nuc subpencil Auc-λEuc contains the finite eigenvalues of M - λN (also called the uncontrollable finite eigenvalues of A - λE). If E is singular, Auc-λEuc may also contain a part of the infinite eigenvalues of M - λN (also called the uncontrollable infinite eigenvalues of A - λE).

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 B, and the relative tolerance for the nonzero elements of A and B. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary left transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed orthogonal or unitary right transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

Method: The implemented algorithm [1] represents a specialization of the controllability staircase algorithm of [2] to the special structure of the input matrix B = [B1,B2].

References

[1] Varga, A. Reliable algorithms for computing minimal dynamic covers for descriptor systems. Proc. of MTNS'04, Leuven, Belgium, 2004.

[2] Varga, A. Computation of Irreducible Generalized State-Space Realizations. Kybernetika, vol. 26, pp. 89-106, 1990.

source
MatrixPencils.sklf_leftfin!Function
sklf_leftfin!(A, E, C, B; fast = true, atol1 = 0, atol2 = 0,  
              rtol, withQ = true, withZ = true) -> (Q, Z, μl, no, nuo)

Reduce the partitioned full column rank matrix pencil

  M - λN = | A-λE | n
           |  C   | p
              n

with A-λE regular, to an equivalent form F - λG = Q'*(M - λN)*diag(I,Z) using orthogonal or unitary transformation matrices Q and Z such that M - λN is transformed into a staircase form exhibiting the separation of its finite eigenvalues:

                | Auo-λEuo  |   *    | nuo
  | At-λEt |    |   0       | Ao-λEo | no
  |--------| =  |-----------|--------|
  |  Ct    |    |   0       |   Co   | p
                    nuo         no

Ct = C*Z, At = Q'*A*Z and Et = Q'*E*Z are returned in C, A and E, respectively, and Q'*B is returned in B (unless B = missing). The resulting Et is upper triangular. If E is already upper triangular, then the preliminary reduction of E to upper triangular form is not performed.

The subpencil | Ao-λEo | has full column rank no for all finite values of λ, is in a staircase form, | Co | and, if E is nonsingular, contains the left Kronecker indices of M - λN. The nl-dimensional vector μl contains the row and column dimensions of the blocks of the staircase form such that i-th block has dimensions μl[nl-i] x μl[nl-i+1] (with μl[0] = p) and has full column rank. If E is nonsingular, the difference μl[nl-i]-μl[nl-i+1] for i = 1, 2, ..., nl is the number of elementary Kronecker blocks of size i x (i-1).

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 C, and the relative tolerance for the nonzero elements of A and C. The reduction is performed using rank decisions based on rank revealing QR-decompositions with column pivoting if fast = true or the more reliable SVD-decompositions if fast = false.

The performed orthogonal or unitary left transformations are accumulated in the matrix Q if withQ = true. Otherwise, Q is set to nothing. The performed orthogonal or unitary right transformations are accumulated in the matrix Z if withZ = true. Otherwise, Z is set to nothing.

Note: This function, called with reversed input parameters E and A (i.e., instead A and E), performs the separation all infinite and nonzero finite eigenvalues of the pencil M - λN.

source