Manipulation of general linear matrix pencils

FunctionDescription
pbalance!Balancing arbitrary matrix pencils.
pbalqualBalancing quality of a matrix pencils.
klfComputation of the Kronecker-like form exhibiting the full Kronecker structure
klf_rightComputation of the Kronecker-like form exhibiting the right and finite Kronecker structures
klf_rightinfComputation of the Kronecker-like form exhibiting the right and infinite Kronecker structures
klf_leftComputation of the Kronecker-like form exhibiting the left and finite Kronecker structures
klf_leftinfComputation of the Kronecker-like form exhibiting the left and infinite Kronecker structures
klf_rlsplitComputation of the Kronecker-like form exhibiting the separation of right and left Kronecker structures
preduceBFReduction to the basic condensed form [B A-λE; D C] with E upper triangular and nonsingular.
klf_right!Computation of the Kronecker-like form exhibiting the right and finite Kronecker structures
klf_right_refine!Update the Kronecker-like form by splitting the right and infinite Kronecker structures
klf_right_refineut!Refine the Kronecker-like form by enforcing upper triangular shapes of blocks in the leading full row rank subpencil
klf_right_refineinf!Refine the Kronecker-like form by enforcing upper triangular shapes of blocks in its infinite regular part
klf_left!Computation of the Kronecker-like form exhibiting the left and finite Kronecker structures
klf_left_refine!Update the Kronecker-like form by splitting the left and infinite Kronecker structures
klf_left_refineut!Refine the Kronecker-like form by enforcing upper triangular shapes of blocks in the leading full row rank subpencil
klf_left_refineinf!Refine the Kronecker-like form by enforcing upper triangular shapes of blocks in its infinite regular part
MatrixPencils.pbalance!Function
pbalance!(M, N; r, c, regpar, shift, maxiter, tol, pow2) -> (Dl, Dr)

Balance the m×n matrix pencil M - λN by reducing the 1-norm of the matrix T := abs(M)+abs(N) by row and column balancing. This involves similarity transformations with diagonal matrices Dl and Dr applied to T to make the rows and columns of Dl*T*Dr as close in norm as possible. The modified Sinkhorn–Knopp algorithm described in [1] is employed to reduce T to an approximately doubly stochastic matrix. The targeted row and column sums can be specified using the keyword arguments r = rs and c = cs, where rs and cs are m- and n-dimensional positive vectors, representing the desired row and column sums, respectively (Default: rs = ones(m) and cs = ones(n)).

The resulting Dl and Dr are diagonal scaling matrices. If the keyword argument pow2 = true is specified, then the components of the resulting optimal Dl and Dr are replaced by their nearest integer powers of 2. If pow2 = false, the optimal values Dl and Dr are returned. The resulting Dl*M*Dr and Dl*N*Dr overwrite M and N, respectively

A regularization-based scaling is performed if a nonzero regularization parameter α is specified via the keyword argument regpar = α. If diagreg = true, then the balancing algorithm is performed on the extended symmetric matrix [ α^2*I T; T' α^2*I ], while if diagreg = false (default), the balancing algorithm is performed on the matrix [ (α/m)^2*em*em' T; T' (α/n)^2*en*en' ], where em and en are m- and n-dimensional vectors with elements equal to one. If α = 0 and shift = γ > 0 is specified, then the algorithm is performed on the rank-one perturbation T+γ*em*en.

The keyword argument tol = τ, with τ ≤ 1, specifies the tolerance used in the stopping criterion. The iterative process is stopped as soon as the incremental scalings are tol-close to the identity.

The keyword argument maxiter = k specifies the maximum number of iterations k allowed in the balancing algorithm.

Method: This function employs the regularization approaches proposed in [1], modified to handle matrices with zero rows or zero columns. The alternative shift based regularization has been proposed in [2].

[1] F.M.Dopico, M.C.Quintana and P. van Dooren, "Diagonal scalings for the eigenstructure of arbitrary pencils", SIMAX, 43:1213-1237, 2022.

[2] P.A.Knight, The Sinkhorn–Knopp algorithm: Convergence and applications, SIAM J. Matrix Anal. Appl., 30 (2008), pp. 261–275.

source
MatrixPencils.pbalqualFunction
qs = pbalqual(M, N)

Compute the 1-norm based scaling quality of a matrix pencil M-λN.

The resulting qs is computed as

    qs = qS(abs(M)+abs(N)) ,

where qS(⋅) is the scaling quality measure defined in Definition 5.5 of [1] for nonnegative matrices. This definition has been extended to also cover matrices with zero rows or columns. If N = I, qs = qs(M) is computed.

A large value of qs indicates a possibly poorly scaled matrix pencil.

[1] F.M.Dopico, M.C.Quintana and P. van Dooren, "Diagonal scalings for the eigenstructure of arbitrary pencils", SIMAX, 43:1213-1237, 2022.

source
MatrixPencils.klfFunction
klf(M, N; 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 matrix pencil M - λN 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 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.klf_rightFunction
klf_right(M, N; fast = true, ut = false, atol1 = 0, atol2 = 0, rtol, withQ = true, withZ = true) -> (F, G, Q, Z, νr, μr, nf, ν, μ)

Reduce the matrix pencil M - λN 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.klf_rightinfFunction
klf_rightinf(M, N; fast = true, ut = false, atol1 = 0, atol2 = 0, rtol, withQ = true, withZ = true) 
            -> (F, G, Q, Z, νr, μr, νi, n, p)

Reduce the matrix pencil M - λN 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 infinite and left Kronecker structures:

               |  Mr-λNr  |     *      |    *     |
               |----------|------------|----------|
     F - λG =  |    O     |   Mi-λNi   |    *     |
               |----------|------------|----------|
               |    O     |     0      | Mfl-λNfl |

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 regular pencil Mi-λNi is in a staircase form, contains the infinite elementary divisors of M-λN with Mi upper triangular if ut = true and nonsingular, and Ni is upper triangular and nilpotent. The ni-dimensional vector νi contains the dimensions of the square diagonal blocks of the staircase form Mi-λNi such that the i-th block has dimensions νi[i] x νi[i]. 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 Mfl-λNfl contains the left Kronecker and finite Kronecker structures of the pencil M-λN and is in the form

                 | A-λE | 
      Mfl-λNfl = |------| ,        
                 |  C   |

where E is an nxn non-singular upper triangular matrix, and A and C are nxn- and pxn-dimensional matrices, respectively.

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.klf_rlsplitFunction
klf_rlsplit(M, N; fast = true, finite_infinite = false, atol1 = 0, atol2 = 0, rtol, withQ = true, withZ = true) -> (F, G, Q, Z, ν, μ, n, m, p)

Reduce the matrix pencil M - λN 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 one of the following Kronecker-like forms:

(1) if finite_infinite = false, then

               | Mri-λNri |     *      |
      F - λG = |----------|------------|, 
               |    O     | Mfl-λNfl   |

where the subpencil Mri-λNri contains the right Kronecker structure and infinite elementary divisors and the subpencil Mfl-λNfl contains the finite and left Kronecker structure of the pencil M-λN.

The full row rank subpencil Mri-λNri is in a staircase form. 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).

The full column rank subpencil Mfl-λNfl is in the form

                 | A-λE | 
      Mfl-λNfl = |------| ,        
                 |  C   |

where E is an nxn non-singular upper triangular matrix, and A and C are nxn- and pxn-dimensional matrices, respectively, and m = 0.

(2) if finite_infinite = true, then

               | Mrf-λNrf |     *      |
      F - λG = |----------|------------|, 
               |    O     | Mil-λNil   |

where the subpencil Mrf-λNrf contains the right Kronecker and finite Kronecker structures and the subpencil Mil-λNil contains the left Kronecker structures and infinite elementary divisors of the pencil M-λN.

The full row rank subpencil Mrf-λNrf is in the form

      Mrf-λNrf = | B  | A-λE | ,

where E is an nxn non-singular upper triangular matrix, and A and B are nxn- and nxm-dimensional matrices, respectively, and p = 0.

The full column rank sub pencil Mil-λNil is in a staircase form. 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).

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.preduceBFFunction
preduceBF(M, N; fast = true, atol = 0, rtol, withQ = true, withZ = true) -> (F, G, Q, Z, n, m, p)

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

               | B  | A-λE | 
      F - λG = |----|------| ,        
               | D  |  C   |

where E is an nxn non-singular matrix, and A, B, C, D are nxn-, nxm-, pxn- and pxm-dimensional matrices, respectively. The order n of E is equal to the numerical rank of N determined using the absolute tolerance atol and relative tolerance rtol. M and N are overwritten by F and G, respectively.

If fast = true, E is determined upper triangular using a rank revealing QR-decomposition with column pivoting of N 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, E is determined diagonal using a rank revealing SVD-decomposition of N 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.klf_right!Function
klf_right!(M, N; fast = true, roff = 0, coff = 0, rtrail = 0, ctrail = 0, atol = 0, rtol, withQ = true, withZ = true) -> (Q, Z, νr, μr, nf, ν, μ, tol)

Reduce the partitioned matrix pencil M - λN (* stands for a not relevant subpencil)

          [   *         *        *  ] roff
 M - λN = [   0      M22-λN22    *  ] m
          [   0         0        *  ] rtrail
            coff        n     ctrail

to an equivalent form F - λG = Q1'*(M - λN)*Z1 using orthogonal or unitary transformation matrices Q1 and Z1 such that the subpencil M22 - λN22 is transformed into the following Kronecker-like form exhibiting its right and finite Kronecker structures:

               |  Mr-λNr    |     *      |    *     |
               |------------|------------|----------|
 F22 - λG22 =  |    O       |   Mf-λNf   |    *     |
               |------------|------------|----------|
               |    O       |     0      | Mil-λNil |

F and G are returned in M and N, respectively.

The full row rank pencil Mr-λNr, in a staircase form, contains the right Kronecker indices of the subpencil M22 - λN22 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.

The nf x nf pencil Mf-λNf is regular and contains the finite elementary divisors of the subpencil M22 - λN22. 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 subpencil M22 - λN22. 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 ν[nb-i+1] x μ[nb-i+1] and has full column rank. The difference ν[nb-i+1]-μ[nb-i+1] for i = 1, 2, ..., nb is the number of elementary Kronecker blocks of size i x (i-1). The difference μ[nb-i+1]-ν[nb-i] for i = 1, 2, ..., nb is the number of infinite elementary divisors of degree i (with ν[0] = 0).

The keyword arguments atol and rtol, specify the absolute and relative tolerances for the nonzero elements of M, respectively. The internally employed absolute tolerance for the nonzero elements of M is returned in tol. 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 (i.e., Q <- Q*Q1) if withQ = true. Otherwise, Q is not modified. The performed right orthogonal or unitary transformations are accumulated in the matrix Z (i.e., Z <- Z*Z1) if withZ = true. Otherwise, Z is not modified.

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

source
MatrixPencils.klf_right_refine!Function
 klf_right_refine!(ν, μ, M, N, tol; fast = true, ut = false, roff = 0, coff = 0, rtrail = 0, ctrail = 0, withQ = true, withZ = true) -> (νr, μr, νi)

Reduce the partitioned matrix pencil M - λN (* stands for a not relevant subpencil)

         [   *         *        *  ] roff
M - λN = [   0      Mri-λNri    *  ] mri
         [   0         0        *  ] rtrail
           coff      nri     ctrail

to an equivalent form F - λG = Q1'*(M - λN)*Z1 using orthogonal or unitary transformation matrices Q1 and Z1 such that the full row rank subpencil Mri-λNriis transformed into the following Kronecker-like form exhibiting its right and infinite Kronecker structures:

              |  Mr-λNr    |     *      |
 Fri - λGri = |------------|------------|
              |    O       |   Mi-λNi   |

The full row rank pencil Mri-λNri is in a staircase form and the nb-dimensional vectors ν and μ contain the row and, respectively, column dimensions of the blocks of the staircase form Mri-λNri such that the i-th block has dimensions ν[i] x μ[i] and has full row rank. The matrix Mri has full row rank and the trailing μ[2]+...+μ[nb] columns of Nri form a full column rank submatrix.

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] = μ[i]-ν[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 regular pencil Mi-λNiis in a staircase form, contains the infinite elementary divisors of Mri-λNri, Mi is upper triangular if ut = true and nonsingular and Ni is upper triangular and nilpotent. The ni-dimensional vector νi contains the dimensions of the square diagonal blocks of the staircase form Mi-λNi such that the i-th block has dimensions νi[i] x νi[i]. The difference νi[ni-i+1]-νi[ni-i] = ν[i]-μ[i+1] for i = 1, 2, ..., ni is the number of infinite elementary divisors of degree i (with νi[0] = 0 and μ[nb+1] = 0).

F and G are returned in M and N, respectively.

The performed left orthogonal or unitary transformations are accumulated in the matrix Q (i.e., Q <- Q*Q1) if withQ = true. Otherwise, Q is not modified. The performed right orthogonal or unitary transformations are accumulated in the matrix Z (i.e., Z <- Z*Z1) if withZ = true. Otherwise, Z is not modified.

source
MatrixPencils.klf_right_refineut!Function
klf_right_refineut!(ν, μ, M, N, Q, Z; ctrail = 0, withQ = true, withZ = true)

Reduce the partitioned matrix pencil M - λN (* stands for a not relevant subpencil)

         [ Mri-λNri  *   ] 
M - λN = [    0      *   ] 
                   ctrail

to an equivalent form F - λG = Q1'*(M - λN)*Z1 using orthogonal or unitary transformation matrices Q1 and Z1, such that the full row rank subpencil Mri-λNri, in staircase form , is transformed as follows: 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. It is assumed that Mri-λNri has nb diagonal blocks and the dimensions of the diagonal blocks are specified by the nb-dimensional vectors ν and μ such that the i-th block has dimensions ν[i] x μ[i].

The performed orthogonal or unitary transformations are accumulated in Q (i.e., Q <- Q*Q1), if withQ = true and Z (i.e., Z <- Z*Z1), if withZ = true.

source
MatrixPencils.klf_right_refineinf!Function
klf_right_refineinf!(νi, M, N, Z, R; roff = 0, coff = 0, withZ = true)

Reduce the partitioned matrix pencil M - λN (* stands for a not relevant subpencil)

         [   *         *        *  ] roff
M - λN = [   0       Mi-λNi     *  ] ni
         [   0         0        *  ] rtrail
           coff       ni     ctrail

to an equivalent form F - λG = (M - λN)*Z1 using an orthogonal or unitary transformation matrix Z1, such that the regular subpencil Mi-λNi, in staircase form , with Mi nonsingular and Ni nillpotent and upper triangular, is transformed to obtain Mi is upper-triangular and preserve Ni upper triangular. It is assumed that Mi-λNi has nb diagonal blocks and the i-th diagonal block has dimensions νi[i] x νi[i].

The performed orthogonal or unitary transformations are accumulated in Z (i.e., Z <- Z*Z1), if withZ = true.

The matrix R is overwritten by R*Z1 unless R = missing.

source
MatrixPencils.klf_leftFunction
klf_left(M, N; fast = true, ut = false, atol1 = 0, atol2 = 0, rtol, withQ = true, withZ = true) -> (F, G, Q, Z, ν, μ, nf, νl, μl)

Reduce the matrix pencil M - λN 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). 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.

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.klf_leftinfFunction
klf_leftinf(M, N; fast = true, ut = false, atol1 = 0, atol2 = 0, rtol, withQ = true, withZ = true) 
            -> (F, G, Q, Z, n, m, νi, νl, μl)

Reduce the matrix pencil M - λN 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 infinite and left Kronecker structures:

               |  Mrf-λNrf  |     *      |    *    |
               |------------|------------|---------|
     F - λG =  |    O       |   Mi-λNi   |    *    |
               |------------|------------|---------|
               |    O       |     0      |  Ml-λNl |

The full row rank pencil Mrf-λNrf contains the right Kronecker and finite Kronecker structures of the pencil M-λN and is in the standard form

      Mrf-λNrf = | B  | A-λE | ,

where E is an nxn non-singular upper triangular matrix, and A and B are nxn- and nxm-dimensional matrices, respectively.

The regular pencil Mi-λNi is in a staircase form, contains the infinite elementary divisors of M-λN with Mi upper triangular if ut = true and nonsingular, and Ni is upper triangular and nilpotent. The ni-dimensional vector νi contains the dimensions of the square diagonal blocks of the staircase form Mi-λNi such that the i-th block has dimensions νi[i] x νi[i]. The difference νi[i]-νi[i-1] = ν[ni-i+1]-μ[ni-i+1] for i = 1, 2, ..., ni is the number of infinite elementary divisors of degree i (with νi[0] = 0 and μ[nb+1] = 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.klf_left!Function
klf_left!(M, N; fast = true, roff = 0, coff = 0, rtrail = 0, ctrail = 0, atol = 0, rtol, withQ = true, withZ = true) -> (Q, Z, ν, μ, nf, νl, μl, tol)

Reduce the partitioned matrix pencil M - λN (* stands for a not relevant subpencil)

          [   *         *        *  ] roff
 M - λN = [   0      M22-λN22    *  ] m
          [   0         0        *  ] rtrail
            coff        n     ctrail

to an equivalent form F - λG = Q1'*(M - λN)*Z1 using orthogonal or unitary transformation matrices Q1 and Z1 such that the subpencil M22 - λN22 is transformed into the following Kronecker-like form exhibiting its finite and left Kronecker structures

               |  Mri-λNri  |     *      |    *    |
               |------------|------------|---------|
 F22 - λG22 =  |    O       |   Mf-λNf   |    *    |
               |------------|------------|---------|
               |    O       |     0      |  Ml-λNl |

F and G are returned in M and N, respectively.

The full row rank pencil Mri-λNri is in a staircase form, and contains the right Kronecker indices and infinite elementary divisors of the subpencil M22 - λN22. 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).

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 nl-dimensional vectors νl and μl contain the row and, respectively, column dimensions of the blocks of the staircase form Ml-λNl such that j-th block has dimensions νl[nl-j+1] x μl[nl-j+1] 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 atol and rtol, specify the absolute and relative tolerances for the nonzero elements of M, respectively. The internally employed absolute tolerance for the nonzero elements of M is returned in tol. 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 (i.e., Q <- Q*Q1) if withQ = true. Otherwise, Q is not modified. The performed right orthogonal or unitary transformations are accumulated in the matrix Z (i.e., Z <- Z*Z1) if withZ = true. Otherwise, Z is not modified.

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

source
MatrixPencils.klf_left_refine!Function
klf_left_refine!(ν, μ, M, N, tol; fast = true, ut = false, roff = 0, coff = 0, rtrail = 0, ctrail = 0, withQ = true, withZ = true) -> (νi, νl, μl)

Reduce the partitioned matrix pencil M - λN (* stands for a not relevant subpencil)

         [   *         *        *  ] roff
M - λN = [   0      Mil-λNil    *  ] mli
         [   0         0        *  ] rtrail
           coff      nli     ctrail

to an equivalent form F - λG = Q1'*(M - λN)*Z1 using orthogonal or unitary transformation matrices Q1 and Z1 such that the full column rank subpencil Mil-λNil is transformed into the following Kronecker-like form exhibiting its infinite and left Kronecker structures:

              |  Mi-λNi    |     *      |
 Fil - λGil = |------------|------------|
              |    O       |   Ml-λNl   |

The full column rank pencil Mil-λNil is in a staircase form and 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 matrix Mil has full column rank and the leading μ[1]+...+μ[nb-1] columns of Nil form a full row rank submatrix.

The regular pencil Mi-λNi is in a staircase form, contains the infinite elementary divisors of Mil-λNil, Mi is upper triangular if ut = true and nonsingular and Ni is upper triangular and nilpotent. The nbi-dimensional vector νi contains the dimensions of the square diagonal blocks of the staircase form Mi-λNi such that the i-th block has dimensions νi[i] x νi[i]. The difference νi[i]-νi[i-1] = ν[nbi-i+1]-μ[nbi-i+1] for i = 1, 2, ..., nbi is the number of infinite elementary divisors of degree i (with νi[0] = 0 and μ[nbi+1] = 0).

The full column rank pencil Mil-λNil is in a staircase form, and contains the left Kronecker indices and infinite elementary divisors of the subpencil M22 - λN22. 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 ν[nb-i+1] x μ[nb-i+1] and has full column rank. The difference ν[nb-i+1]-μ[nb-i+1] for i = 1, 2, ..., nb is the number of elementary Kronecker blocks of size i x (i-1). The difference μ[nb-i+1]-ν[nb-i] for i = 1, 2, ..., nb is the number of infinite elementary divisors of degree i (with ν[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] = ν[nl-j+1]-μ[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 [Y; 0] with Y 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.

F and G are returned in M and N, respectively.

The performed left orthogonal or unitary transformations are accumulated in the matrix Q (i.e., Q <- Q*Q1) if withQ = true. Otherwise, Q is not modified. The performed right orthogonal or unitary transformations are accumulated in the matrix Z (i.e., Z <- Z*Z1) if withZ = true. Otherwise, Z is not modified.

source
MatrixPencils.klf_left_refineut!Function
klf_left_refineut!(ν, μ, M, N, Q, Z; roff = 0, coff = 0, withQ = true, withZ = true)

Reduce the partitioned matrix pencil M - λN (* stands for a not relevant subpencil)

         [  *     *     ] roff
M - λN = [  0  Mil-λNil ] 
          coff

to an equivalent form F - λG = Q1'*(M - λN)*Z1 using orthogonal or unitary transformation matrices Q1 and Z1, such that the full column rank subpencil Mil-λNil, in staircase form , is transformed as follows: the full column rank diagonal blocks of Mil are reduced to the form [Y ; 0] with Y 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. It is assumed that Mil-λNil has nb diagonal blocks and the dimensions of the diagonal blocks are specified by the nb-dimensional vectors ν and μ such that the i-th block has dimensions ν[i] x μ[i].

The performed orthogonal or unitary transformations are accumulated in Q (i.e., Q <- Q*Q1), if withQ = true and Z (i.e., Z <- Z*Z1), if withZ = true.

source
MatrixPencils.klf_left_refineinf!Function
klf_left_refineinf!(νi, M, N, Q, L; roff = 0, coff = 0, withQ = true)

Reduce the partitioned matrix pencil M - λN (* stands for a not relevant subpencil)

         [   *         *        *  ] roff
M - λN = [   0       Mi-λNi     *  ] ni
         [   0         0        *  ] rtrail
           coff       ni     ctrail

to an equivalent form F - λG = Q1'*(M - λN) using an orthogonal or unitary transformation matrix Q1, such that the regular subpencil Mi-λNi, in staircase form, with Mi nonsingular and Ni nillpotent and upper triangular, is transformed to obtain Mi upper-triangular and preserve Ni upper triangular. It is assumed that Mi-λNi has nb diagonal blocks and the i-th diagonal block has dimensions νi[i] x νi[i].

The performed orthogonal or unitary transformations are accumulated in Q (i.e., Q <- Q*Q1), if withQ = true.

The matrix L is overwritten by Q1'*L unless L = missing.

source