Manipulation of polynomial matrices

FunctionDescription
poly2pmConversion of a polynomial matrix used in Polynomials package to a polynomial matrix represented as a 3-dimensional matrix
pm2polyConversion of a polynomial matrix represented as a 3-dimensional matrix to a polynomial matrix used in Polynomials package
pmdegDetermination of the degree of a polynomial matrix
pmevalEvaluation of a polynomial matrix for a given value of its argument.
pmreverseBuilding the reversal of a polynomial matrix
pmdivremQuotients and remainders of elementwise divisions of two polynomial matrices
pm2lpCF1Building a linearization in the first companion Frobenius form
pm2lpCF2Building a linearization in the second companion Frobenius form
pm2lsBuilding a descriptor system based structured linearization [A-λE B; C D] of a polynomial matrix
ls2pmComputation of the polynomial matrix from its descriptor system based structured linearization
pm2lpsBuilding a pencil based structured linearization [A-λE B-λF; C-λG D-λH] of a polynomial matrix
lps2pmComputation of the polynomial matrix from its pencil based structured linearization
spm2lsBuilding a descriptor system based structured linearization [A-λE B; C D] of a structured polynomial matrix [T(λ) U(λ); V(λ) W(λ)]
spm2lpsBuilding a pencil based structured linearization [A-λE B-λF; C-λG D-λH] of a structured polynomial matrix [T(λ) U(λ); V(λ) W(λ)]
MatrixPencils.poly2pmFunction
poly2pm(PM; grade = k) -> P

Build a grade k matrix polynomial representation P(λ) from a polynomial matrix, polynomial vector or scalar polynomial PM(λ).

PM(λ) is a matrix, vector or scalar of elements of the Polynomial type provided by the Polynomials package.

P(λ) is a grade k polynomial matrix of the form P(λ) = P_1 + λ P_2 + ... + λ**k P_(k+1), for which the coefficient matrices P_l, l = 1, ..., k+1, are stored in the 3-dimensional matrix P, where P[:,:,l] contains the l-th coefficient matrix P_l (multiplying λ**(l-1)). If grade = missing, then k is chosen the largest degree of the elements of PM. The coefficients of the degree d element (i,j) of PM(λ) result in P[i,j,1:d+1].

source
MatrixPencils.pm2polyFunction
pm2poly(P[,var = 'x']) -> PM

Build the polynomial matrix PM(λ) from its matrix polynomial representation P(λ).

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

PM(λ) is a matrix of elements of the Polynomial type provided by the Polynomials package. The element (i,j) of PM(λ) is built from the coefficients contained in P[i,j,1:k+1]. The symbol to be used for the indeterminate λ can be specified in the optional input variable var.

source
MatrixPencils.pmdegFunction
pmdeg(P) -> deg

Determine the degree deg of a polynomial matrix P(λ).

P(λ) is a grade k polynomial matrix of the form P(λ) = P_1 + λ P_2 + ... + λ**k P_(k+1), for which the coefficient matrices P_i, i = 1, ..., k+1, are stored in the 3-dimensional matrix P, where P[:,:,i] contains the i-th coefficient matrix P_i (multiplying λ**(i-1)). The degree of P(λ) is deg = j-1, where j is the largest index for which P[:,:,j] is nonzero. The degree of the zero polynomial matrix is defined to be deg = -1.

P(λ) can also be specified as a matrix, vector or scalar of elements of the Polynomial type provided by the Polynomials package. The degree of P(λ) is the largest degree of the elements of P(λ). The degree of the zero polynomial matrix is defined to be -1.

source
MatrixPencils.pmevalFunction
pmeval(P,val) -> R

Evaluate R = P(val) for a polynomial matrix P(λ), using Horner's scheme.

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

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

source
MatrixPencils.pmreverseFunction
pmreverse(P[,j]) -> Q

Build Q(λ) = λ^j*P(1/λ), the j-reversal of a polynomial matrix P(λ) for j ≥ deg(P(λ)). If j is not specified, the default value j = deg(P(λ)) is used.

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

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

If deg(P(λ)), then Q(λ) is a grade j polynomial matrix of the form Q(λ) = Q_1 + λ Q_2 + ... + λ**j Q_(j+1), for which the coefficient matrices Q_i, i = 1, ..., j+1, are stored in the 3-dimensional matrix Q, where Q[:,:,i] contains the i-th coefficient matrix Q_i (multiplying λ**(i-1)). The coefficient matrix Q_i is either 0 if i ≤ j-d or Q_(j-d+i) = P_(d-i+1) for i = 1, ..., d.

source
MatrixPencils.pmdivremFunction
pmdivrem(N,D) -> (Q, R)

Compute the quotients in Q(λ) and remainders in R(λ) of the elementwise polynomial divisions N(λ)./D(λ).

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

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

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

The polynomial matrices of quotients Q(λ) and remainders R(λ) are stored in the 3-dimensional matrices Q and R, respectively, where Q[:,:,i] and R[:,:,i] contain the i-th coefficient matrix multiplying λ**(i-1).

source
MatrixPencils.pm2lpCF1Function
 pm2lpCF1(P; grade = l) -> (M, N)

Build a strong linearization M - λN of a polynomial matrix P(λ) in the first companion Frobenius form.

P(λ) is a grade k polynomial matrix assumed of the form P(λ) = P_1 + λ P_2 + ... + λ**k P_(k+1), with the coefficient matrices P_i, i = 1, ..., k+1 stored in the 3-dimensional matrix P, where P[:,:,i] contains the i-th coefficient matrix P_i (multiplying λ**(i-1)). The effective grade l to be used for linearization can be specified via the keyword argument grade as grade = l, where l must be chosen equal to or greater than the degree of P(λ). The default value used for l is l = deg(P(λ)).

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

If P(λ) is a m x n polynomial matrix of effective grade l and degree d, then the resulting matrix pencil M - λN satisfies the following conditions [1]:

(1) M - λN has dimension (m+n*(l-1)) x n*l and M - λN is regular if P(λ) is regular;

(2) M - λN and P(λ) have the same finite eigenvalues;

(3) the partial multiplicities of infinite eigenvalues of M - λN are in excess with l-d to the partial multiplicities of the infinite eigenvalues of P(λ);

(4) M - λN and P(λ) have the same number of right Kronecker indices and the right Kronecker indices of M - λN are in excess with l-1 to the right Kronecker indices of P(λ);

(5) M - λN and P(λ) have the same left Kronecker structure (i.e., the same left Kronecker indices).

[1] F. De Terán, F. M. Dopico, D. S. Mackey, Spectral equivalence of polynomial matrices and the Index Sum Theorem, Linear Algebra and Its Applications, vol. 459, pp. 264-333, 2014.

source
MatrixPencils.pm2lpCF2Function
pm2lpCF2(P; grade = l) -> (M, N)

Build a strong linearization M - λN of a polynomial matrix P(λ) in the second companion Frobenius form.

P(λ) is a grade k polynomial matrix assumed of the form P(λ) = P_1 + λ P_2 + ... + λ**k P_(k+1), with the coefficient matrices P_i, i = 1, ..., k+1 stored in the 3-dimensional matrix P, where P[:,:,i] contains the i-th coefficient matrix P_i (multiplying λ**(i-1)). The effective grade l to be used for linearization can be specified via the keyword argument grade as grade = l, where l must be chosen equal to or greater than the degree of P(λ). The default value used for l is l = deg(P(λ)).

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

If P(λ) is a m x n polynomial matrix of effective grade l and degree d, then the resulting matrix pencil M - λN satisfies the following conditions [1]:

(1) M - λN has dimension l*m x (n+(l-1)*m) and M - λN is regular if P(λ) is regular;

(2) M - λN and P(λ) have the same finite eigenvalues;

(3) the partial multiplicities of infinite eigenvalues of M - λN are in excess with l-d to the partial multiplicities of the infinite eigenvalues of P(λ);

(4) M - λN and P(λ) have the same right Kronecker structure (i.e., the same right Kronecker indices);

(5) M - λN and P(λ) have the same number of left Kronecker indices and the left Kronecker indices of M - λN are in excess with l-1 to the left Kronecker indices of P(λ).

[1] F. De Terán, F. M. Dopico, D. S. Mackey, Spectral equivalence of polynomial matrices and the Index Sum Theorem, Linear Algebra and Its Applications, vol. 459, pp. 264-333, 2014.

source
MatrixPencils.pm2lsFunction
 pm2ls(P; contr = false, obs = false, noseig = false, minimal = false,
          fast = true, atol = 0, rtol) -> (A, E, B, C, D)

Build a structured linearization as a system matrix S(λ) of the form

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

of the polynomial matrix P(λ) which preserves a part of the Kronecker structure of P(λ).

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

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

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

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

(2) rank[B A-λE] = n (controllability) if minimal = true or contr = true, in which case the right Kronecker structure is preserved;

(3) rank[A-λE; C] = n (observability) if minimal = true or contr = true, in which case the left Kronecker structure is preserved;

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

If conditions (1)-(4) are satisfied, the linearization is called minimal and the resulting order n is the least achievable order. If conditions (1)-(3) are satisfied, the linearization is called irreducible and the resulting order n is the least achievable order using orthogonal similarity transformations. For an irreducible linearization S(λ) preserves the pole-zero structure (finite and infinite) and the left and right Kronecker structures of P(λ).

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

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

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

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

source
MatrixPencils.pm2lpsFunction
 pm2lps(P; contr = false, obs = false) -> (A, E, B, F, C, G, D, H)

Build a structured linearization as a system matrix S(λ) of the form

        | A-λE | B-λF | 
 S(λ) = |------|------|
        | C-λG | D-λH |

of the polynomial matrix P(λ) which preserves a part of the Kronecker structure of P(λ).

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

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

If d is the degree of the p x m polynomial matrix P(λ), then the computed linearization satisfies:

(1) A-λE is a n x n regular pencil, where n = p(d-1) if contr = false and p <= m and n = m(d-1) otherwise;

(2) P(λ) = (C-λG)*inv(λE-A)*(B-λF)+D-λH;

(3) rank[B-λF A-λE] = n for any finite and infinite λ (strong controllability) if contr = true, in which case the right Kronecker structure is preserved;

(4) rank[A-λE; C-λG] = n for any finite and infinite λ (strong observability) if obs = true, in which case the left Kronecker structure is preserved.

If conditions (1)-(4) are satisfied, the linearization is called strongly minimal, the resulting order n is the least achievable order and S(λ) preserves the pole-zero structure (finite and infinite) and the left and right Kronecker structures of P(λ).

The pencil based linearization is built using the methods described in [1].

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

source
MatrixPencils.spm2lsFunction
 spm2ls(T, U, V, W; fast = true, contr = false, obs = false, minimal = false, atol = 0, rtol) -> (A, E, B, C, D)

Build a structured linearization as a system matrix S(λ) of the form

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

of the structured polynomial matrix

          | -T(λ) | U(λ) |
   P(λ) = |-------|------|
          | V(λ)  | W(λ) |

such that V(λ)*inv(T(λ))*U(λ)+W(λ) = C*inv(λE-A)*B+D. The resulting linearization S(λ) preserves a part, if minimal = false, or the complete Kronecker structure, if minimal = true, of P(λ). In the latter case, the order n of A-λE is the least possible one and S(λ) is a strong linearization of P(λ).

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

T(λ), U(λ), V(λ), and W(λ) can also be specified as matrices, vectors or scalars of elements of the Polynomial type provided by the Polynomials package.

The computed structured linearization satisfies:

(1) A-λE is regular;

(2) rank[B A-λE] = n (controllability) if minimal = true or contr = true, in which case the finite and right Kronecker structures are preserved;

(3) rank[A-λE; C] = n (observability) if minimal = true or obs = true, in which case the finite and left Kronecker structures are preserved;

(4) A-λE has no simple infinite eigenvalues if minimal = true, in which case the complete Kronecker structure is preserved.

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

The structured linearization is built using the methods described in [1].

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

source
MatrixPencils.spm2lpsFunction
 spm2lps(T, U, V, W; fast = true, contr = false, obs = false, minimal = false, atol = 0, rtol) -> (A, E, B, F, C, G, D, H)

Build a structured linearization

          | A-λE | B-λF | 
 M - λN = |------|------|
          | C-λG | D-λH |

of the structured polynomial matrix

          | -T(λ) | U(λ) |
   P(λ) = |-------|------|
          |  V(λ) | W(λ) |

such that V(λ)*inv(T(λ))*U(λ)+W(λ) = (C-λG))*inv(λE-A)*(B-λF)+D-λH. The resulting linearization M - λN preserves a part, if minimal = false, or the complete Kronecker structure, if minimal = true, of P(λ). In the latter case, the order n of A-λE is the least possible one and M - λN is a strong linearization of P(λ).

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

T(λ), U(λ), V(λ), and W(λ) can also be specified as matrices, vectors or scalars of elements of the Polynomial type provided by the Polynomials package.

The computed structured linearization satisfies:

(1) A-λE is regular;

(2) rank[B-λF A-λE] = n (strong controllability) if minimal = true or contr = true, in which case the finite and right Kronecker structures are preserved;

(3) rank[A-λE; C-λG] = n (strong observability) if minimal = true or obs = true, in which case the finite and left Kronecker structures are preserved.

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

source
MatrixPencils.ls2pmFunction
ls2pm(A, E, B, C, D; fast = true, atol1 = 0, atol2 = 0, gaintol = 0, rtol = min(atol1,atol2) > 0 ? 0 : n*ϵ, val) -> P

Build the polynomial matrix P(λ) = C*inv(λE-A)*B+D corresponding to its structured linearization

 | A-λE | B | 
 |------|---|
 |  C   | D |

by explicitly determining for each polynomial entry, its coefficients from its roots and a corresponding gain.

The keyword arguments atol1 and atol2 specify the absolute tolerances for the elements of A, B, C, D, and, respectively, of E, and rtol specifies the relative tolerances for the nonzero elements of A, B, C, D and E. The default relative tolerance is (n+1)*ϵ, where n is the size of the size dimension of A, and ϵ is the machine epsilon of the element type of coefficients of A.

The keyword argument gaintol specifies the threshold for the magnitude of the nonzero elements of the gain matrix C*inv(γE-A)*B+D, where γ = val if val is a number or γ is a randomly chosen complex value of unit magnitude, if val = missing. Generally, val should not be a root of any of entries of P.

source
MatrixPencils.lps2pmFunction
lps2pm(A, E, B, F, C, G, D, H; fast = true, atol1 = 0, atol2 = 0, gaintol = 0, rtol = min(atol1,atol2) > 0 ? 0 : n*ϵ, val) -> P

Build the polynomial matrix P(λ) = (C-λG)*inv(λE-A)*(B-λF)+D-λH corresponding to its structured linearization

 | A-λE | B-λF | 
 |------|------|
 | C-λG | D-λH |

by explicitly determining for each polynomial entry, its coefficients from its roots and corresponding gain.

The keyword arguments atol1 and atol2 specify the absolute tolerances for the elements of A, B, C, D, and of E, F, G, H, respectively, and rtol specifies the relative tolerances for the nonzero elements of A, B, C, D, E, F,G,H. The default relative tolerance is(n+2)*ϵ, wherenis the size of the size dimension ofA, andϵis the machine epsilon of the element type of coefficients ofA`.

The keyword argument gaintol specifies the threshold for the magnitude of the nonzero elements of the gain matrix C*inv(γE-A)*B+D, where γ = val if val is a number or γ is a randomly chosen complex value of unit magnitude, if val = missing. Generally, val should not be a root of any of entries of P.

source