Some applications of matrix pencil computations
Function | Description |
---|---|
pkstruct | Determination of the complete Kronecker structure |
prank | Determination of the normal rank |
peigvals | Computation of the finite and infinite eigenvalues |
pzeros | Computation of the finite and infinite zeros |
MatrixPencils.KRInfo
— TypeKRInfo
Kronecker-structure object definition.
If info::KRInfo
is the Kronecker-structure object, then:
info.rki
is a vector whose components contains the right Kronecker indices;
info.lki
is a vector whose components contains the left Kronecker indices;
info.id
is a vector whose components contains the orders of the infinite elementary divisors (i.e., the multiplicities of infinite eigenvalues).
info.nf
is the number of finite eigenvalues.
info.nrank
is the normal rank.
Destructuring via iteration produces the components info.rki
, info.lki
, info.id
, info.nf
and info.nrank
.
MatrixPencils.pkstruct
— Functionpkstruct(M, N; fast = false, atol1::Real = 0, atol2::Real = 0, rtol::Real=min(atol1,atol2)>0 ? 0 : n*ϵ) -> KRInfo
Determine the Kronecker-structure information of the matrix pencil M-λN
and return an KRInfo
object.
The right Kronecker indices rki
, left Kronecker indices lki
, infinite elementary divisors id
, the number of finite eigenvalues nf
and the normal rank nrank
can be obtained from KRInfo
as KRInfo.rki
, KRInfo.lki
, KRInfo.id
, KRInfo.nf
and KRInfo.nrank
, respectively. The determination of the Kronecker-structure information is performed by reducing the pencil M-λN
to an appropriate Kronecker-like form (KLF) exhibiting all structural elements of the pencil M-λN
. The reduction is performed using orthogonal similarity transformations and involves rank decisions based on rank revealing QR-decompositions with column pivoting, if fast = true
, or, the more reliable, SVD-decompositions, if fast = false
. For efficiency purposes, the reduction is only partially performed, without accumulating the performed orthogonal transformations.
The right Kronecker indices are provided in the integer vector rki
, where each i
-th element rki[i]
is the row dimension k
of an elementary Kronecker block of size k x (k+1)
. The number of elements of rki
is the dimension of the right nullspace of the pencil M-λN
and their sum is the degree of a right minimal polynomial nullspace basis.
The left Kronecker indices are provided in the integer vector lki
, where each i
-th element lki[i]
is the column dimension k
of an elementary Kronecker block of size (k+1) x k
. The number of elements of lki
is the dimension of the left nullspace of the pencil M-λN
and their sum is the degree of a left minimal polynomial nullspace basis.
The multiplicities of infinite eigenvalues are provided in the integer vector id
, where each i
-th element id[i]
is the order of an infinite elementary divisor (i.e., the multiplicity of an infinite eigenvalue). The number of elements of id
is the number of infinite elementary divisors of the pencil M-λN
.
The normal rank nrank
is the maximum rank of M-λN
for all values of λ
. If N = nothing
, then nrank
is the rank of M
.
The keyword arguements atol1
, atol2
and rtol
specify 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
, respectively. The default relative tolerance is n*ϵ
, where n
is the size of the smallest dimension of M
, and ϵ
is the machine epsilon of the element type of M
.
MatrixPencils.prank
— Functionprank(M::AbstractMatrix, N::AbstractMatrix; fastrank = true, atol1::Real=0, atol2::Real=0, rtol::Real=min(atol1,atol2)>0 ? 0 : n*ϵ)
Compute the normal rank of a matrix pencil M - λN
. If fastrank = true
, the rank is evaluated by counting how many singular values of M - γ N
have magnitude greater than max(max(atol1,atol2), rtol*σ₁)
, where σ₁
is the largest singular value of M - γ N
and γ
is a randomly generated value. If fastrank = false
, the rank is evaluated as nr + ni + nf + nl
, where nr
and nl
are the sums of right and left Kronecker indices, respectively, while ni
and nf
are the number of infinite and finite eigenvalues, respectively. The sums nr+ni
and nf+nl
, are determined from an appropriate Kronecker-like form (KLF) exhibiting the spliting of the right and left structures of the pencil M - λN
.
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 default relative tolerance is n*ϵ
, where n
is the size of the smallest dimension of M
, and ϵ
is the machine epsilon of the element type of M
. For efficiency purpose, the reduction to the relevant KLF is only partially performed using rank decisions based on rank revealing SVD-decompositions.
The use of atol
and rtol
keyword arguments in rank determinations requires at least Julia 1.1. To enforce compatibility with Julia 1.0, the newer function rank in Julia 1.1 has been explicitly included.
MatrixPencils.peigvals
— Functionpeigvals(M, N; fast = false, atol1::Real = 0, atol2::Real = 0, rtol::Real=min(atol1,atol2)>0 ? 0 : n*ϵ) -> (val, KRInfo)
Return the (finite and infinite) eigenvalues of the matrix pencil M-λN
in val
and information on the Kronecker-structure in the KRInfo
object.
The information on the Kronecker-structure consists of the right Kronecker indices rki
, left Kronecker indices lki
, infinite elementary divisors id
, the number of finite eigenvalues nf
, and the normal rank nrank
and can be obtained from KRInfo
as KRInfo.rki
, KRInfo.lki
, KRInfo.id
, KRInfo.nf
and KRInfo.nrank
, respectively. For more details, see pkstruct
.
The computation of the eigenvalues is performed by reducing the pencil M-λN
to an appropriate Kronecker-like form (KLF) exhibiting the spliting of the infinite and finite eigenvalue structures of the pencil M-λN
. The reduction is performed using orthonal similarity transformations and involves rank decisions based on rank revealing QR-decompositions with column pivoting, if fast = true
, or, the more reliable, SVD-decompositions, if fast = false
. For efficiency purposes, the reduction is only partially performed, without accumulating the performed orthogonal transformations.
The keyword arguements atol1
, atol2
and rtol
specify 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
, respectively. The default relative tolerance is n*ϵ
, where n
is the size of the smallest dimension of M
, and ϵ
is the machine epsilon of the element type of M
.
MatrixPencils.pzeros
— Functionpzeros(M, N; fast = false, atol1::Real = 0, atol2::Real = 0, rtol::Real=min(atol1,atol2)>0 ? 0 : n*ϵ) -> (val, iz, KRInfo)
Return the finite and infinite zeros of the matrix pencil M-λN
in val
, the multiplicities of infinite zeros in iz
and information on the Kronecker-structure in the KRInfo
object.
The information on the multiplicities of infinite zeros is provided in the vector iz
, where each i
-th element iz[i]
is equal to k-1
, where k
is the order of an infinite elementary divisor with k > 0
. The number of infinite zeros contained in val
is the sum of the components of iz
.
The information on the Kronecker-structure consists of the right Kronecker indices rki
, left Kronecker indices lki
, infinite elementary divisors id
, the number of finite eigenvalues nf
, and the normal rank nrank
and can be obtained from KRInfo
as KRInfo.rki
, KRInfo.lki
, KRInfo.id
, KRInfo.nf
and KRInfo.nrank
, respectively. For more details, see pkstruct
.
The computation of the zeros is performed by reducing the pencil M-λN
to an appropriate Kronecker-like form (KLF) exhibiting the spliting of the infinite and finite eigenvalue structures of the pencil M-λN
. The reduction is performed using orthonal similarity transformations and involves rank decisions based on rank revealing QR-decompositions with column pivoting, if fast = true
, or, the more reliable, SVD-decompositions, if fast = false
. For efficiency purposes, the reduction is only partially performed, without accumulating the performed orthogonal transformations.
The keyword arguements atol1
, atol2
and rtol
specify 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
, respectively. The default relative tolerance is n*ϵ
, where n
is the size of the smallest dimension of M
, and ϵ
is the machine epsilon of the element type of M
.