Periodic Schur decompositions
phess
Periodic Hessenberg decomposition of a product of matrices.phess!
Periodic Hessenberg decomposition of a product of matrices (in place computation).pschur
Periodic Schur decompositions of products or quotient products of matrices.pschur!
Periodic Schur decompositions of products of matrices (in place computation).psordschur!
Reordering of periodic Schur decompositions of products or quotient products of matrices.psordschur1!
Reordering of periodic Schur decompositions of products or quotient products of square matrices.pgschur
Generalized real periodic Schur decomposition of a formal product of matrices.pgschur!
Generalized real periodic Schur decompositions of formal products of matrices (in place computation).pgordschur!
Reordering of generalized real periodic Schur decompositions a formal products of matrices.
PeriodicMatrices.phess
— Function phess(A::Array{Float64,3}; hind = 1, rev = true, withZ = true) -> (H, Z, ihess)
phess1(A::Array{Float64,3}; hind = 1, rev = true, withZ = true) -> (H, Z, ihess)
Compute the Hessenberg decomposition of a product of square matrices A(p)*...*A(2)*A(1)
, if rev = true
(default) or A(1)*A(2)*...*A(p)
if rev = false
, without evaluating the product. The matrices A(1)
, ...
, A(p)
are contained in the n×n×p
array A
such that the i
-th matrix A(i)
is contained in A[:,:,i]
. The resulting n×n×p
arrays H
and Z
contain the matrices H(1)
, ...
, H(p)
and the orthogonal matrices Z(1)
, ...
, Z(p)
, respectively, such that for rev = true
Z(2)' * A(1) * Z(1) = H(1),
Z(3)' * A(2) * Z(2) = H(2),
...
Z(1)' * A(p) * Z(p) = H(p),
and for rev = false
Z(1)' * A(1) * Z(2) = H(1),
Z(2)' * A(2) * Z(3) = H(2),
...
Z(p)' * A(p) * Z(1) = H(p).
If hind = ihess
, with 1 ≤ ihess ≤ p
(default ihess = 1
), then H(i)
, i = 1, ..., p
are in a periodic Hessenberg form, with H(ihess)
in upper Hessenberg form and H(i)
upper triangular for i
$\neq$ ihess
. H(i)
and Z(i)
are contained in H[:,:,i]
and Z[:,:,i]
, respectively. The performed orthogonal transformations are not accumulated if withZ = false
, in which case Z = nothing
.
The function phess
is based on a wrapper for the SLICOT subroutine MB03VW
(see PeriodicMatrices.SLICOTtools.mb03vw!
).
The function phess1
is based on wrappers for the SLICOT subroutines MB03VD
(see PeriodicMatrices.SLICOTtools.mb03vd!
) and MB03VY
(see PeriodicMatrices.SLICOTtools.mb03vy!
).
PeriodicMatrices.phess!
— Function phess!(A::Array{Float64,3}, ilh::Tuple(Int,Int) = (1, size(A,1)); kwargs...) -> (H, Z, ihess)
Same as phess(A; kwargs...)
but uses the input matrix A
as workspace and specifies a range ilh = (ilo, ihi)
, such that all matrices A(j), j = 1, ..., p
, are already in periodic Hessenberg forms in rows and columns 1:ilo-1
and ihi+1:n
, where n
is the first dimension of A
.
PeriodicMatrices.pschur
— Function pschur(A::Array{Float64,3}; rev = true, withZ = true) -> (S, Z, ev, ischur, α, γ)
pschur1(A::Array{Float64,3}; rev = true, withZ = true) -> (S, Z, ev, ischur, α, γ)
pschur2(A::Array{Float64,3}; rev = true, withZ = true) -> (S, Z, ev, ischur, α, γ)
Compute the Schur decomposition of a product of square matrices A(p)*...*A(2)*A(1)
, if rev = true
(default) or A(1)*A(2)*...*A(p)
if rev = false
, without evaluating the product. The matrices A(1)
, ...
, A(p)
are contained in the n×n×p
array A
such that the i
-th matrix A(i)
is contained in A[:,:,i]
. The resulting n×n×p
arrays S
and Z
contain the matrices S(1)
, ...
, S(p)
and the orthogonal matrices Z(1)
, ...
, Z(p)
, respectively, such that for rev = true
Z(2)' * A(1) * Z(1) = S(1),
Z(3)' * A(2) * Z(2) = S(2),
...
Z(1)' * A(p) * Z(p) = S(p),
and for rev = false
Z(1)' * A(1) * Z(2) = S(1),
Z(2)' * A(2) * Z(3) = S(2),
...
Z(p)' * A(p) * Z(1) = S(p).
If sind = ischur
, with 1 ≤ ischur ≤ p
(default ischur = 1
), then S(i)
, for i = 1, ..., p
are in a periodic Schur form, with S(ischur)
in quasi-upper triangular (or Schur) form and S(i)
upper triangular for i
$\neq$ ischur
. S(i)
and Z(i)
are contained in S[:,:,i]
and Z[:,:,i]
, respectively. The vector ev
contains the eigenvalues of the appropriate matrix product. The eigenvalues can be alternatively expressed as α .* γ
, where γ
contains suitable scaling parameters to avoid overflows or underflows in the expressions of the eigenvalues. The performed orthogonal transformations are not accumulated if withZ = false
, in which case Z = nothing
.
The function pschur
is based on wrappers for the SLICOT subroutines MB03VW
(see PeriodicMatrices.SLICOTtools.mb03vw!
) and MB03BD
(see PeriodicMatrices.SLICOTtools.mb03bd!
), based on algorithms proposed in [1] and [2].
The function pschur1
is based on wrappers for the SLICOT subroutines MB03VD
(see PeriodicMatrices.SLICOTtools.mb03vd!
), MB03VY
(see PeriodicMatrices.SLICOTtools.mb03vy!
) and MB03BD
(see PeriodicMatrices.SLICOTtools.mb03bd!
), based on algorithms proposed in [1] and [2].
The function pschur2
is based on wrappers for the SLICOT subroutines MB03VD
(see PeriodicMatrices.SLICOTtools.mb03vd!
), MB03VY
(see PeriodicMatrices.SLICOTtools.mb03vy!
) and MB03VW
(see PeriodicMatrices.SLICOTtools.mb03vw!
), based on the algorithm proposed in [1]. Known issue: MB03VW
may fails for larger periods.
REFERENCES
[1] Bojanczyk, A., Golub, G. H. and Van Dooren, P. The periodic Schur decomposition: algorithms and applications. In F.T. Luk (editor), Advanced Signal Processing Algorithms, Architectures, and Implementations III, Proc. SPIE Conference, vol. 1770, pp. 31-42, 1992.
[2] Kressner, D. An efficient and reliable implementation of the periodic QZ algorithm. IFAC Workshop on Periodic Control Systems (PSYCO 2001), Como (Italy), August 27-28 2001. Periodic Control Systems 2001 (IFAC Proceedings Volumes), Pergamon.
pschur(A::Vector{Matrix{T}}; rev = true, withZ = true) -> (S, Z, ev, ischur, α, γ)
pschur1(A::Vector{Matrix{T}}; rev = true, withZ = true) -> (S, Z, ev, ischur, α, γ)
Compute the extended periodic Schur decomposition of a square product of matrices A(p)*...*A(2)*A(1)
, if rev = true
(default) or A(1)*A(2)*...*A(p)
if rev = false
, without evaluating the product. The matrices A(1)
, ...
, A(p)
are contained in the p
-vector of matrices A
such that the i
-th matrix A(i)
, of dimensions m(i)×n(i)
, is contained in A[i]
. The resulting p
-vectors S
and Z
contain the matrices S(1)
, ...
, S(p)
and the orthogonal matrices Z(1)
, ...
, Z(p)
, respectively, such that for rev = true
Z(2)' * A(1) * Z(1) = S(1),
Z(3)' * A(2) * Z(2) = S(2),
...
Z(1)' * A(p) * Z(p) = S(p),
and for rev = false
Z(1)' * A(1) * Z(2) = S(1),
Z(2)' * A(2) * Z(3) = S(2),
...
Z(p)' * A(p) * Z(1) = S(p).
The resulting index ischur
is determined such that m(ischur) ≤ m(i), ∀i
. The resulting S(i)
, for i = 1, ..., p
are in an extended periodic Schur form, with S(ischur)
in a quasi-upper trapezoidal form and S(i)
upper trapezoidal for i
$\neq$ ischur
. S(i)
and Z(i)
are contained in S[i]
and Z[i]
, respectively. The first nmin
components of ev := α .* γ
contain the core eigenvalues of the appropriate matrix product, where nmin = m(ischur)
, while the last nmax-nmin
components of ev
are zero, where nmax
is the largest row or column dimension of A(i)
, for i = 1, ..., p
. The eigenvalues can be alternatively expressed as α .* γ
, where γ
contains suitable scaling parameters to avoid overflows or underflows in the expressions of the eigenvalues. The performed orthogonal transformations are not accumulated if withZ = false
, in which case Z = nothing
.
The function pschur
is based on wrappers for the SLICOT subroutines MB03VW
(see PeriodicMatrices.SLICOTtools.mb03vw!
) and MB03BD
(see PeriodicMatrices.SLICOTtools.mb03bd!
), based on algorithms proposed in [1] and [2].
The function pschur1
is based on wrappers for the SLICOT subroutines MB03VD
(see PeriodicMatrices.SLICOTtools.mb03vd!
), MB03VY
(see PeriodicMatrices.SLICOTtools.mb03vy!
), and MB03BD
(see PeriodicMatrices.SLICOTtools.mb03bd!
), based on algorithms proposed in [1] and [2].
REFERENCES
[1] Bojanczyk, A., Golub, G. H. and Van Dooren, P. The periodic Schur decomposition: algorithms and applications. In F.T. Luk (editor), Advanced Signal Processing Algorithms, Architectures, and Implementations III, Proc. SPIE Conference, vol. 1770, pp. 31-42, 1992.
[2] Kressner, D. An efficient and reliable implementation of the periodic QZ algorithm. IFAC Workshop on Periodic Control Systems (PSYCO 2001), Como (Italy), August 27-28 2001. Periodic Control Systems 2001 (IFAC Proceedings Volumes), Pergamon.
pschur(A::AbstractArray{T,3}, E::AbstractArray{T,3}; rev = true, withZ = true) -> (S, T, Q, Z, ev, ischur, α, γ)
Compute the periodic Schur decomposition of a square formal quotient product of matrices inv(E(p))*A(p)*...*inv(E(2))*A(2)*inv(E(1))*A(1)
, if rev = true
(default) or A(1)*inv(E(1))*A(2)*inv(E(2))*...*A(p)*inv(E(p))
if rev = false
, without evaluating the product and the inverses. The matrices A(1)
, ...
, A(p)
are contained in the n×n×p
-array A
such that the i
-th matrix A(i)
is contained in A[:,:,i]
. The square matrices E(1)
, ...
, E(p)
are contained in the n×n×p
-array E
such that the i
-th matrix E(i)
is contained in E[:,:,i]
.
The resulting n×n×p
-arrays S
, T
, Q
and Z
contain, respectively, the matrices S(1)
, ...
, S(p)
with S(ischur)
in a quasi-upper trapezoidal form and S(i)
upper trapezoidal for i
$\neq$ ischur
, the upper triangular matrices T(1)
, ...
, T(p)
, the orthogonal matrices Q(1)
, ...
, Q(p)
, and Z(1)
, ...
, Z(p)
, such that for rev = true
Q(1)' * A(1) * Z(1) = S(1), Q(1)' * E(1) * Z(2) = T(1),
Q(2)' * A(2) * Z(2) = S(2), Q(2)' * E(2) * Z(3) = T(2),
...
Q(p)' * A(p) * Z(p) = S(p), Q(p)' * E(p) * Z(1) = T(p),
and for rev = false
Q(1)' * A(1) * Z(1) = S(1), Q(2)' * E(1) * Z(1) = T(1),
Q(2)' * A(2) * Z(2) = S(2), Q(3)' * E(2) * Z(2) = T(2),
...
Q(p)' * A(p) * Z(p) = S(p), Q(1)' * E(p) * Z(p) = T(p),
The complex vector ev
contains the eigenvalues of the appropriate matrix product, and can be alternatively expressed as ev := α .* γ
, where γ
contains suitable scaling parameters to avoid overflows or underflows in the expressions of the eigenvalues. The performed orthogonal transformations are not accumulated if withZ = false
, in which case Q = nothing
and Z = nothing
.
The function pschur
is based on wrappers for the SLICOT subroutines MB03VW
(see PeriodicMatrices.SLICOTtools.mb03vw!
) and MB03BD
(see PeriodicMatrices.SLICOTtools.mb03bd!
), based on algorithms proposed in [1] and [2].
REFERENCES
[1] Bojanczyk, A., Golub, G. H. and Van Dooren, P. The periodic Schur decomposition: algorithms and applications. In F.T. Luk (editor), Advanced Signal Processing Algorithms, Architectures, and Implementations III, Proc. SPIE Conference, vol. 1770, pp. 31-42, 1992.
[2] Kressner, D. An efficient and reliable implementation of the periodic QZ algorithm. IFAC Workshop on Periodic Control Systems (PSYCO 2001), Como (Italy), August 27-28 2001. Periodic Control Systems 2001 (IFAC Proceedings Volumes), Pergamon.
pschur(A::Vector{Matrix{T}}, E::Vector{Matrix{T}}; rev = true, withZ = true) -> (S, T, Q, Z, ev, ischur, α, γ)
Compute the extended periodic Schur decomposition of a square formal product of matrices inv(E(p))*A(p)*...*inv(E(2))*A(2)*inv(E(1))*A(1)
, if rev = true
(default) or A(1)*inv(E(1))*A(2)*inv(E(2))*...*A(p)*inv(E(p))
if rev = false
, without evaluating the product and the inverses. The matrices A(1)
, ...
, A(p)
are contained in the p
-vector of matrices A
such that the i
-th matrix A(i)
, of dimensions m(i)×n(i)
, is contained in A[i]
. The square matrices E(1)
, ...
, E(p)
are contained in the p
-vector of matrices E
such that the i
-th matrix E(i)
, of dimensions m(i)×m(i)
if rev = true
or n(i)×n(i)
if rev = false
, is contained in E[i]
.
The resulting index ischur
is determined such that m(ischur) ≤ m(i), ∀i
. The resulting p
-vectors S
, T
, Q
and Z
contain, respectively, the matrices S(1)
, ...
, S(p)
with S(ischur)
in a quasi-upper trapezoidal form and S(i)
upper trapezoidal for i
$\neq$ ischur
, the upper triangular matrices T(1)
, ...
, T(p)
, the orthogonal matrices Q(1)
, ...
, Q(p)
, and Z(1)
, ...
, Z(p)
, such that for rev = true
Q(1)' * A(1) * Z(1) = S(1), Q(1)' * E(1) * Z(2) = T(1),
Q(2)' * A(2) * Z(2) = S(2), Q(2)' * E(2) * Z(3) = T(2),
...
Q(p)' * A(p) * Z(p) = S(p), Q(p)' * E(p) * Z(1) = T(p),
and for rev = false
Q(1)' * A(1) * Z(1) = S(1), Q(2)' * E(1) * Z(1) = T(1),
Q(2)' * A(2) * Z(2) = S(2), Q(3)' * E(2) * Z(2) = T(2),
...
Q(p)' * A(p) * Z(p) = S(p), Q(1)' * E(p) * Z(p) = T(p),
The first nmin
components of ev := α .* γ
contain the core eigenvalues of the appropriate matrix product, where nmin = m(ischur)
, while the last nmax-nmin
components of ev
are zero, where nmax
is the largest row or column dimension of A(i)
, for i = 1, ..., p
. The eigenvalues can be alternatively expressed as α .* γ
, where γ
contains suitable scaling parameters to avoid overflows or underflows in the expressions of the eigenvalues. The performed orthogonal transformations are not accumulated if withZ = false
, in which case Q = nothing
and Z = nothing
.
The function pschur
is based on wrappers for the SLICOT subroutines MB03VW
(see PeriodicMatrices.SLICOTtools.mb03vw!
) and MB03BD
(see PeriodicMatrices.SLICOTtools.mb03bd!
), based on algorithms proposed in [1] and [2].
REFERENCES
[1] Bojanczyk, A., Golub, G. H. and Van Dooren, P. The periodic Schur decomposition: algorithms and applications. In F.T. Luk (editor), Advanced Signal Processing Algorithms, Architectures, and Implementations III, Proc. SPIE Conference, vol. 1770, pp. 31-42, 1992.
[2] Kressner, D. An efficient and reliable implementation of the periodic QZ algorithm. IFAC Workshop on Periodic Control Systems (PSYCO 2001), Como (Italy), August 27-28 2001. Periodic Control Systems 2001 (IFAC Proceedings Volumes), Pergamon.
PeriodicMatrices.pschur!
— Function pschur!(A::Array{Float64,3}, ilh::Tuple(Int,Int) = (1, size(A,1)); kwargs...) -> (S, Z, ihess)
Same as pschur(A; kwargs...)
but uses the input matrix A
as workspace and specifies a range ilh = (ilo, ihi)
, such that all matrices A(j), j = 1, ..., p
, are already in periodic Schur forms in rows and columns 1:ilo-1
and ihi+1:n
, where n
is the first dimension of A
.
pschur!(A::Array{Float64,3}, Z::AbstractArray{Float64,3}, ilh::Tuple(Int,Int) = (1, size(A,1)); rev = true, sind = 1, withZ = true) -> (ev, sind, α, γ)
pschur!(wspace::Tuple, A::Array{Float64,3}, Z::AbstractArray{Float64,3}, ilh::Tuple(Int,Int) = (1, size(A,1)); rev = true, sind = 1, withZ = true) -> (ev, sind, α, γ)
Compute the Schur decomposition of a product of square matrices A(p)*...*A(2)*A(1)
, if rev = true
(default) or A(1)*A(2)*...*A(p)
if rev = false
, without evaluating the product. The matrices A(1)
, ...
, A(p)
are contained in the n×n×p
array A
such that the i
-th matrix A(i)
is contained in A[:,:,i]
. A range ilh = (ilo, ihi)
can be specified, such that all matrices A(j), j = 1, ..., p
, are already in periodic Schur forms in rows and columns 1:ilo-1
and ihi+1:n
, where n
is the first dimension of A
. The resulting reduced matrices S(1)
, ...
, S(p)
, representing the periodic Schur decimposition overwrite the input matrices A(1)
, ...
, A(p)
. If withZ = true
, then Z
must contain on input an n×n×p
array, which will be overwritten on output with the orthogonal matrices Z(1)
, ...
, Z(p)
used for reduction. If withZ = false
, orthogonal transformations are not computed and Z
can be an arbiitrary three-dimensional array (e.g., an empty 0×0×0
array). The resulting reduced matrices S(1)
, ...
, S(p)
and the orthogonal transformation matrices Z(1)
, ...
, Z(p)
, satisfy for rev = true
Z(2)' * A(1) * Z(1) = S(1),
Z(3)' * A(2) * Z(2) = S(2),
...
Z(1)' * A(p) * Z(p) = S(p),
and for rev = false
Z(1)' * A(1) * Z(2) = S(1),
Z(2)' * A(2) * Z(3) = S(2),
...
Z(p)' * A(p) * Z(1) = S(p).
If sind = ischur
, with 1 ≤ ischur ≤ p
(default ischur = 1
), then S(i)
, for i = 1, ..., p
are in a periodic Schur form, with S(ischur)
in quasi-upper triangular (or Schur) form and S(i)
upper triangular for i
$\neq$ ischur
. The vector ev
contains the eigenvalues of the appropriate matrix product. The eigenvalues can be alternatively expressed as α .* γ
, where γ
contains suitable scaling parameters to avoid overflows or underflows in the expressions of the eigenvalues.
To reduce the number of allocations when the periodic Schur decomposition needs to be repeteadly computed, the internally needed workspace can be pre-allocated in a tuple wspace
using the function PeriodicMatrices.ws_pschur
).
PeriodicMatrices.ws_pschur
— Function ws_pschur(A::Array{Float64,3}) -> (QIND, SIND, ALPHAR, ALPHAI, BETA, SCAL, IWORK, DWORK)
Allocate workspaces for the function pschur!
to call the SLICOT-based wrappers mb03vw!
and mb03bd!
avoiding additional allocations.
PeriodicMatrices.psordschur!
— Function psordschur!(S::Vector{Matrix{Float64}}, Z::Vector{Matrix{Float64}}, select; rev, schurindex)
Reorder the core eigenvalues of the product Π = S[p]*...*S[2]*S[1]
, if rev = true
(default) or Π = S[1]*S[2]*...*S[p]
if rev = false
, where Π
is in real Schur form, such that the selected eigenvalues in the logical array select
are moved into the leading positions. The p
-vectors S
and Z
contain the matrices S[1]
, ...
, S[p]
in an extended periodic Schur form, with the leading square block of S[schurindex]
in real Schur form, and the corresponding orthogonal transformation matrices Z[1]
, ...
, Z[p]
, respectively. S
and Z
are overwritten by the updated matrices. A conjugate pair of eigenvalues must be either both included or both excluded via select
. The dimension of select must be equal to the number of core eigenvalues (i.e., the minimum dimension of matrices in the vector S
).
psordschur!(S::Array{Float64,3}, Z::Array{Float64,3}, select::Union{Vector{Bool},BitVector}; rev, schurindex)
Reorder the eigenvalues of the product Π = S[:,:,p]*...*S[:,:,2]*S[:,:,1]
, if rev = true
(default) or Π = S[:,:,1]*S[:,:,2]*...*S[:,:,p]
if rev = false
, where Π
is in real Schur form, such that the selected eigenvalues in the logical array select
are moved into the leading positions. The 3-dimensional arrays S
and Z
contain the matrices S[:,:,1]
, ...
, S[:,:,p]
in a periodic Schur form, with S[:,:,schurindex]
in real Schur form, and the corresponding orthogonal transformation matrices Z[:,:,1]
, ...
, Z[:,:,p]
, respectively. S
and Z
are overwritten by the updated matrices. A conjugate pair of eigenvalues must be either both included or both excluded via select
.
PeriodicMatrices.psordschur1!
— Function psordschur1!(S::Vector{Matrix{Float64}}, Z::Vector{Matrix{Float64}}, select; rev, schurindex)
Reorder the eigenvalues of the product Π = S[p]*...*S[2]*S[1]
, if rev = true
(default) or Π = S[1]*S[2]*...*S[p]
if rev = false
, where Π
is in real Schur form, such that the selected eigenvalues in the logical array select
are moved into the leading positions. The p
-vectors S
and Z
contain, respectively, the square matrices S[1]
, ...
, S[p]
in a periodic Schur form, with S[schurindex]
in real Schur form, and the corresponding orthogonal transformation matrices Z[1]
, ...
, Z[p]
, respectively. S
and Z
are overwritten by the updated matrices. A conjugate pair of eigenvalues must be either both included or both excluded via select
.
PeriodicMatrices.pgschur
— Function pgschur(A::Vector{Matrix}, s::Union{Vector{Bool},BitVector}; rev = true, withZ = true) -> (S, Z, ev, ischur, α, γ)
Compute the generalized real periodic Schur decomposition of a formal product of square matrices A(p)^s(p)*...A(2)^s(2)*A(1)^s(1)
, if rev = true
(default), or A(1)^s(1)*A(2)^s(2)*...*A(p)^s(p)
, if rev = false
, where 's(j) = ±1'. The matrices A(1)
, A(2)
, ...
, A(p)
are contained in the p
-dimensional array A
such that the i
-th matrix A(i)
is contained in A[i]
.
The resulting p
-dimensional array S
contains the matrices S(1)
, ...
, S(p)
such that the i
-th matrix S(i)
is contained in S[i]
. The component matrix S[ischur]
is in a quasi-upper triangular form, while S[i]
is upper triangular for i
$\neq$ ischur
. If withZ = true
(default), the resulting p
-dimensional array Z
contains the orthogonal transformation matrices Z(1)
, ...
, Z(p)
such that the i
-th matrix Z(i)
is contained in Z[i]
. The performed orthogonal transformations are not accumulated if withZ = false
, in which case Z = nothing
.
The resulting matrices satisfy for rev = true
Z(mod(j,p)+1)' * A(j) * Z(j) = S(j), if S[j] = true,
Z(j)' * A(j) * Z(mod(j,p)+1) = S(j), if S[j] = false,
and for rev = false
Z(j)' * A(j) * Z(mod(j,p)+1) = S(j), if S[j] = true,
Z(mod(j,p)+1)' * A(j) * Z(j) = S(j), if S[j] = false.
The vector ev
contains the eigenvalues of the appropriate matrix product. The eigenvalues can be alternatively expressed as α .* γ
, where γ
contains suitable scaling parameters to avoid overflows or underflows in the expressions of the eigenvalues.
The function pgschur
is based on wrappers for the SLICOT subroutines MB03VW
(see PeriodicMatrices.SLICOTtools.mb03vw!
) and MB03BD
(see PeriodicMatrices.SLICOTtools.mb03bd!
), based on algorithms proposed in [1] and [2].
REFERENCES
[1] Bojanczyk, A., Golub, G. H. and Van Dooren, P. The periodic Schur decomposition: algorithms and applications. In F.T. Luk (editor), Advanced Signal Processing Algorithms, Architectures, and Implementations III, Proc. SPIE Conference, vol. 1770, pp. 31-42, 1992.
[2] Kressner, D. An efficient and reliable implementation of the periodic QZ algorithm. IFAC Workshop on Periodic Control Systems (PSYCO 2001), Como (Italy), August 27-28 2001. Periodic Control Systems 2001 (IFAC Proceedings Volumes), Pergamon.
pgschur(A::Array{Float64,3}, s::Union{Vector{Bool},BitVector}; rev = true, withQ = true) -> (S, Z, ev, ischur, α, γ)
Compute the generalized real periodic Schur decomposition of a formal product of square matrices A(p)^s(p)*...A(2)^s(2)*A(1)^s(1)
, if rev = true
(default), or A(1)^s(1)*A(2)^s(2)*...*A(p)^s(p)
, if rev = false
, where 's(j) = ±1'. The matrices A(1)
, A(2)
, ...
, A(p)
are contained in the n×n×p
array A
such that the i
-th matrix A(i)
is contained in A[:,:,i]
.
The resulting n×n×p
array S
contains the matrices S(1)
, ...
, S(p)
such that S(ischur)
is in a quasi-upper triangular form, S(i)
is upper triangular for i
$\neq$ ischur
. If withZ = true
(default), the resulting n×n×p
array Z
contains the orthogonal transformation matrices Z(1)
, ...
, Z(p)
. The performed orthogonal transformations are not accumulated if withZ = false
, in which case Z = nothing
.
The resulting matrices satisfy for rev = true
Z(mod(j,p)+1)' * A(j) * Z(j) = S(j), if S[j] = true,
Z(j)' * A(j) * Z(mod(j,p)+1) = S(j), if S[j] = false,
and for rev = false
Z(j)' * A(j) * Z(mod(j,p)+1) = S(j), if S[j] = true,
Z(mod(j,p)+1)' * A(j) * Z(j) = S(j), if S[j] = false.
S(i)
and Z(i)
are contained in S[:,:,i]
and Z[:,:,i]
, respectively. The vector ev
contains the eigenvalues of the appropriate matrix product. The eigenvalues can be alternatively expressed as α .* γ
, where γ
contains suitable scaling parameters to avoid overflows or underflows in the expressions of the eigenvalues.
The function pgschur
is based on wrappers for the SLICOT subroutines MB03VW
(see PeriodicMatrices.SLICOTtools.mb03vw!
) and MB03BD
(see PeriodicMatrices.SLICOTtools.mb03bd!
), based on algorithms proposed in [1] and [2].
REFERENCES
[1] Bojanczyk, A., Golub, G. H. and Van Dooren, P. The periodic Schur decomposition: algorithms and applications. In F.T. Luk (editor), Advanced Signal Processing Algorithms, Architectures, and Implementations III, Proc. SPIE Conference, vol. 1770, pp. 31-42, 1992.
[2] Kressner, D. An efficient and reliable implementation of the periodic QZ algorithm. IFAC Workshop on Periodic Control Systems (PSYCO 2001), Como (Italy), August 27-28 2001. Periodic Control Systems 2001 (IFAC Proceedings Volumes), Pergamon.
PeriodicMatrices.pgschur!
— Functionpgschur!(A::Array{Float64,3}, S::Union{Vector{Bool},BitVector}; kwargs...)
Same as pgschur
but uses the input matrix A
as workspace.
PeriodicMatrices.pgordschur!
— Function pgordschur!(S::Array{Float64,3}, s::Union{Vector{Bool},BitVector}, Z::Array{Float64,3}, select::Union{Vector{Bool},BitVector}; rev, schurindex)
Reorder the eigenvalues of the product Π = S[:,:,p]^s[p]*...*S[:,:,2]^s[2]*S[:,:,1]^s[1]
, if rev = true
(default) or Π = S[:,:,1]^s[1]*S[:,:,2]^s[2]*...*S[:,:,p]^s[p]
if rev = false
, with 's[j] = ±1', where Π
is in a real Schur form, such that the selected eigenvalues in the logical array select
are moved into the leading positions. The 3-dimensional arrays S
and Z
contain the matrices S[:,:,1]
, ...
, S[:,:,p]
in a generalized periodic Schur form, with S[:,:,schurindex]
in a quasi-upper triangular (real Schur) form, and the corresponding orthogonal transformation matrices Z[:,:,1]
, ...
, Z[:,:,p]
, respectively. S
and Z
are overwritten by the updated matrices. A conjugate pair of eigenvalues must be either both included or both excluded via select
.
pgordschur!(S::Vector{Matrix{Float64}}, s::Union{Vector{Bool},BitVector}, Z::Vector{Matrix{Float64}}, select::Union{Vector{Bool},BitVector}; rev, schurindex)
Reorder the eigenvalues of the product Π = S[p]^s[p]*...*S[2]^s[2]*S[1]^s[1]
, if rev = true
(default) or Π = S[1]^s[1]*S[2]^s[2]*...*S[p]^s[p]
if rev = false
, with 's[j] = ±1', where Π
is in a real Schur form, such that the selected eigenvalues in the logical array select
are moved into the leading positions. The p
-vectors S
and Z
contain the matrices S[1]
, ...
, S[p]
in a generalized periodic Schur form, with S[schurindex]
in a quasi-upper triangular (real Schur) form, and the corresponding orthogonal transformation matrices Z[1]
, ...
, Z[p]
, respectively. S
and Z
are overwritten by the updated matrices. A conjugate pair of eigenvalues must be either both included or both excluded via select
.