Constructors for periodic matrices
Discrete-time periodic matrices
PeriodicMatrix
Discrete-time periodic matrix representation.PeriodicArray
Discrete-time periodic array representation.SwitchingPeriodicMatrix
Discrete-time switching periodic matrix representation.SwitchingPeriodicArray
Discrete-time switching periodic array representation.
PeriodicMatrices.PeriodicMatrix
— TypePeriodicMatrix(M, T; nperiod = k) -> A::PeriodicMatrix
Discrete-time periodic matrix representation.
The discrete-time periodic matrix object A
is built from a p
-vector M
of real matrices, the associated time period T
and the number of subperiods specified via the keyword argument nperiod = k
.
M
contains the cyclic component matrices M[i]
, i = 1,..., p
, where M[i]
represents the value M(Δ(i-1))
of a time periodic matrix M(t)
of period T/k
, with Δ := T/(k*p)
, the associated sample time. It is assumed that M[i] := M[mod(i-1,p)+1]
for arbitrary i
. All component matrices are allowed to have arbitrary (time-varying) dimensions. The component matrices M
, the period T
, the number of subperiods k
, the discrete period p
and the sample time Δ
can be accessed via A.M
, A.period
, A.nperiod
, A.dperiod
and A.Ts
, respectively.
PeriodicMatrices.SwitchingPeriodicMatrix
— TypeSwitchingPeriodicMatrix(M, ns, T; nperiod = k) -> A::SwitchingPeriodicMatrix
Discrete-time switching periodic matrix representation.
The discrete-time switching periodic matrix object A
is built from a p
-vector M
of real matrices, a p
-vector ns
of increasing positive integers representing the discrete switching moments, the associated time period T
and the number of subperiods specified via the keyword argument nperiod = k
.
M
contains the component matrices M[i]
, i = 1,..., p
, which defines a sequence of N := ns[p]
of matrices S[1], ..., S[N]
, such that S[j] = M[i]
for j ∈ [ns[i-1]+1, ..., ns[i]]
with ns[0] := 0
. S[j]
is the j
-th value A(Δ(j-1))
of a time periodic matrix A(t)
of subperiod T′ := T/k
, with Δ := T′/N
, the associated sample time. All component matrices are allowed to have arbitrary (time-varying) dimensions. The component matrices M
, the integer vector ns
, the period T
, the number of subperiods k
, the discrete period N
and the sample time Δ
can be accessed via A.M
, A.ns
, A.period
, A.nperiod
, A.dperiod
and A.Ts
, respectively.
The j-th time value A(Δ(j-1))
can be determined as A[j]
. It is assumed that A[j] := A[mod(j-1,N)+1]
for arbitrary j
.
Base.size
— Methodsize(A::PeriodicMatrix)
size(A::SwitchingPeriodicMatrix)
size(A::PeriodicMatrix[, dim])
size(A::SwitchingPeriodicMatrix[, dim])
Return a tuple of two vectors containing the dimensions of the components of the discrete-time periodic matrix A
. Optionally you can specify a dimension dim
to just get the vector of lengths of that dimension.
Base.length
— Methodlength(A::PeriodicMatrix)
length(A::SwitchingPeriodicMatrix)
Return the number of component matrices (also called the discrete period) of the discrete-time periodic matrix A
.
Base.eltype
— Methodeltype(A::PeriodicMatrix)
eltype(A::SwitchingPeriodicMatrix)
Determine the type of the elements of component matrices of the discrete-time periodic matrix A
.
Base.getindex
— Methodgetindex(A::PeriodicMatrix, i)
getindex(A::SwitchingPeriodicMatrix, i)
Return the i
-th component matrix of the discrete-time periodic matrix A
. Equivalent to the syntax A[i]
.
getindex(A::PeriodicMatrix, ind1, ind2)
getindex(A::SwitchingPeriodicMatrix, ind1, ind2)
Return the discrete-time periodic matrix built from the selected ranges [ind1,ind2]
of elements of the component matrices. ind1
and ind2
may be integers, integer ranges or colons.
Base.getindex
— Methodgetindex(A::PeriodicMatrix, ind1, ind2)
getindex(A::SwitchingPeriodicMatrix, ind1, ind2)
Return the discrete-time periodic matrix built from the selected ranges [ind1,ind2]
of elements of the component matrices. ind1
and ind2
may be integers, integer ranges or colons.
Base.lastindex
— Methodlastindex(A::PeriodicMatrix)
lastindex(A::SwitchingPeriodicMatrix)
Return the last index of the component matrices of the discrete-time periodic matrix A
. The syntax A[end]
is equivalent to A[lastindex(A)]
.
Base.lastindex
— Methodlastindex(A::PeriodicMatrix,dim)
lastindex(A::SwitchingPeriodicMatrix,dim)
Return the vector of last indices along dimension dim
of the component matrices of the discrete-time periodic matrix A
.
PeriodicMatrices.PeriodicArray
— TypePeriodicArray(M, T; nperiod = k) -> A::PeriodicArray
Discrete-time periodic array representation.
The discrete-time periodic array object A
is built from a m×n×p
real array M
, the associated time period T
and the number of subperiods specified via the keyword argument nperiod = k
. M
contains the cyclic component matrices M[:,:,i]
, i = 1,..., p
, where M[:,:,i]
represents the value M(Δ(i-1))
of a time periodic matrix M(t)
of period T/k
, with Δ := T/(kp)
, the associated sample time. It is assumed that M[:,:,k] := M[:,:,mod(k-1,p)+1]
for arbitrary k
. The component matrices M
, the period T
, the number of subperiods k
, the discrete period p
and the sample time Δ
can be accessed via A.M
, A.period
, A.nperiod
, A.dperiod
and A.Ts
, respectively.
PeriodicMatrices.SwitchingPeriodicArray
— TypeSwitchingPeriodicArray(M, ns, T; nperiod = k) -> A::SwitchingPeriodicArray
Discrete-time switching periodic array representation.
The discrete-time switching periodic array object A
is built from a m×n×p
real array M
, a p
-vector ns
of increasing positive integers representing the discrete switching moments, the associated time period T
and the number of subperiods specified via the keyword argument nperiod = k
.
M
contains the cyclic component matrices M[:,:,i]
, i = 1,..., p
, which defines a sequence of N := ns[p]
of matrices S[1], ..., S[N]
, such that S[j] =
M[:,:,i]for
j ∈ [ns[i-1]+1, ..., ns[i]]with
ns[0] := 0.
S[j]is the
j-th value
A(Δ(j-1))of a time periodic matrix
A(t)of subperiod
T′ := T/k, with
Δ := T′/N, the associated sample time. The component matrices
M, the integer vector
ns, the period
T, the number of subperiods
k, the discrete period
Nand the sample time
Δcan be accessed via
A.M,
A.ns,
A.period,
A.nperiod,
A.dperiodand
A.Ts`, respectively.
The j-th time value A(Δ(j-1))
can be determined as A[j]
. It is assumed that A[j] := A[mod(j-1,N)+1]
for arbitrary j
.
Base.size
— Methodsize(A::PeriodicArray)
size(A::SwitchingPeriodicArray)
size(A::PeriodicArray[, dim])
size(A::SwitchingPeriodicArray[, dim])
Return a tuple of two integers containing the common row and column dimensions of the components of the discrete-time periodic matrix A
. Optionally you can specify a dimension dim
to just get length of that dimension.
Base.length
— Methodlength(A::PeriodicArray)
length(A::SwitchingPeriodicArray)
Return the number of component matrices (also called the discrete period) of the discrete-time periodic matrix A
.
Base.eltype
— Methodeltype(A::PeriodicArray)
eltype(A::SwitchingPeriodicArray)
Determine the type of the elements of component matrices of the discrete-time periodic matrix A
.
Base.getindex
— Methodgetindex(A::PeriodicArray, i)
getindex(A::SwitchingPeriodicArray, i)
Return the i
-th component matrix of the discrete-time periodic matrix A
. Equivalent to the syntax A[i]
.
Base.getindex
— Methodgetindex(A::PeriodicArray, ind1, ind2)
getindex(A::SwitchingPeriodicArray, ind1, ind2)
Return the discrete-time periodic matrix built from the selected ranges [ind1,ind2]
of elements of the component matrices. ind1
and ind2
may be integers, integer ranges or colons.
Base.lastindex
— Methodlastindex(A::PeriodicArray)
lastindex(A::SwitchingPeriodicArray)
Return the last index of the component matrices of the discrete-time periodic matrix A
. The syntax A[end]
is equivalent to A[lastindex(A)]
.
Base.lastindex
— Methodlastindex(A::PeriodicArray,dim)
lastindex(A::SwitchingPeriodicArray,dim)
Return the last index along dimension dim
of the component matrices of the discrete-time periodic matrix A
.
Continuous-time periodic matrices
PeriodicFunctionMatrix
Continuous-time periodic function matrix representation.PeriodicSymbolicMatrix
Continuous-time periodic symbolic matrix representation.PeriodicTimeSeriesMatrix
Continuous-time periodic time series matrix representation.HarmonicArray
Continuous-time harmonic array representation.FourierFunctionMatrix
Continuous-time Fourier functin matrix representation.PeriodicSwitchingMatrix
Continuous-time switching periodic matrix representation.
PeriodicMatrices.PeriodicFunctionMatrix
— TypePeriodicFunctionMatrix(f, T; nperiod = k) -> A::PeriodicFunctionMatrix
Continuous-time periodic function matrix representation.
The continuous-time periodic function matrix object A
is built from the real matrix function f(t)
of real time variable t
, the associated time period T
and the associated number of subperiods specified via the keyword argument nperiod = k
. It is assumed that f(t) = f(t+T/k)
for any real time value t
. The function f(t)
, the period T
, the row and column dimensions of f(t)
, the number of subperiods k
can be accessed via A.f
, A.period
, the tuple A.dims
and A.nperiod
, respectively.
PeriodicMatrices.PeriodicSymbolicMatrix
— TypePeriodicSymbolicMatrix(F, T; nperiod = k) -> A::PeriodicSymbolicMatrix
Continuous-time periodic symbolic matrix representation.
The continuous-time periodic symbolic matrix object A
is built from F
, a symbolic real matrix or vector of symbolic variable t
, the associated time period T
and the associated number of subperiods specified via the keyword argument nperiod = k
. It is assumed that F(t) = F(t+T/k)
for any real time value t
. The symbolic matrix F
, the period T
and the number of subperiods k
can be accessed via A.F
, A.period
and A.nperiod
, respectively.
PeriodicMatrices.PeriodicTimeSeriesMatrix
— TypePeriodicTimeSeriesMatrix(At, T; nperiod = k) -> A::PeriodicTimeSeriesMatrix
Continuous-time periodic time series matrix representation.
The continuous-time periodic time series matrix object A
of period T
is built from a p
-vector At
of real matrices and the associated subperiod T′ = T/k
, where k ≥ 1
is the number of subperiods (default: k = 1
). At
contains the cyclic component matrices At[i]
, i = 1,..., p
, where At[i]
represents the value A(Δ*(i-1))
of a time periodic matrix A(t)
of period T′
, with Δ := T′/p
, the associated sampling time. It is assumed that At[i] := At[mod(i-1,p)+1]
for arbitrary i
. All component matrices must have the same dimensions. The component matrices At
, the period T
and the number of subperiods k
can be accessed via A.values
, A.period
, and A.nperiod
, respectively.
PeriodicMatrices.HarmonicArray
— Type HarmonicArray(Ahr, T; nperiod = k) -> A::HarmonicArray
Continuous-time harmonic array representation.
The harmonic array object A
of period T
is built using the harmonic representation of a periodic matrix Ahr(t)
of subperiod T′ = T/k
in the form
p
Ahr(t) = A_0 + ∑ ( Ac_i*cos(i*t*2*π/T′)+As_i*sin(i*2*π*t/T′) ) ,
i=1
where k ≥ 1
is the number of subperiods (default: k = 1
). The m×n×(p+1)
complex array Ahr
contains the harmonic components as follows: Ahr[:,:,1]
contains the constant term A_0
(the mean value) and the real and imaginary parts of Ahr[:,:,i+1]
for i = 1, ..., p
contain the coefficient matrices Ac_i
and As_i
, respectively. The complex matrix Ahr
containing the harmonic components, the period T
and the number of subperiods k
can be accessed via A.values
, A.period
and A.nperiod
, respectively.
HarmonicArray(A0, Ac, As, T) -> A::HarmonicArray
Construct a harmonic array representation from the harmonic components.
The harmonic array object A
is built for the harmonic representation Ahr(t)
of a periodic matrix of period T
in the form
p
Ahr(t) = A_0 + ∑ ( Ac_i*cos(i*t*2*π/T)+As_i*sin(i*2*π*t/T) ) ,
i=1
where the constant term A_0
is contained in the real matrix A0
, and Ac
and As
are vectors of real matrices such that the i
-th (cosinus) coefficient matrix Ac_i
is contained in Ac[i]
and the i
-th (sinus) coefficient matrix As_i
is contained in As[i]
. p
is the maximum of length of the vectors of matrices Ac
and As
. If the length of Ac
or As
is less than p
, then zero trailing matrices are assumed in the respective matrix. All component matrices must have the same dimensions. The complex matrix containing the harmonic components and the period T
can be accessed via A.values
and A.period
, respectively.
PeriodicMatrices.HarmonicArray
— Method HarmonicArray(A0, Ac, As, T) -> A::HarmonicArray
Construct a harmonic array representation from the harmonic components.
The harmonic array object A
is built for the harmonic representation Ahr(t)
of a periodic matrix of period T
in the form
p
Ahr(t) = A_0 + ∑ ( Ac_i*cos(i*t*2*π/T)+As_i*sin(i*2*π*t/T) ) ,
i=1
where the constant term A_0
is contained in the real matrix A0
, and Ac
and As
are vectors of real matrices such that the i
-th (cosinus) coefficient matrix Ac_i
is contained in Ac[i]
and the i
-th (sinus) coefficient matrix As_i
is contained in As[i]
. p
is the maximum of length of the vectors of matrices Ac
and As
. If the length of Ac
or As
is less than p
, then zero trailing matrices are assumed in the respective matrix. All component matrices must have the same dimensions. The complex matrix containing the harmonic components and the period T
can be accessed via A.values
and A.period
, respectively.
PeriodicMatrices.FourierFunctionMatrix
— Type FourierFunctionMatrix(Afun, T) -> A::FourierFunctionMatrix
Continuous-time Fourier function matrix representation.
The Fourier function matrix object A
of period T
is built using the Fourier series representation of a periodic matrix Afun(t)
of subperiod T′ = T/k
, where each entry of Afun(t)
has the form
p
a_0 + ∑ ( ac_i*cos(i*t*2*π/T′)+as_i*sin(i*2*π*t/T′) ) ,
i=1
where k ≥ 1
is the number of subperiods (default: k = 1
). The matrix Afun
containing the Fourier representation, the period T
and the number of subperiods k
can be accessed via A.M
, A.period
and A.nperiod
, respectively.
PeriodicMatrices.PeriodicSwitchingMatrix
— TypePeriodicSwitchingMatrix(At, ts, T; nperiod = k) -> A::PeriodicSwitchingMatrix
Continuous-time periodic switching matrix representation.
The continuous-time periodic switching matrix object A
of period T
is built from a p
-vector At
of real matrices, a p
-vector ts
of increasingly ordered switching time values with ts[1] = 0
, and the associated subperiod T′ = T/k
, where k ≥ 1
is the number of subperiods (default: k = 1
). At
contains the cyclic component matrices At[i]
, i = 1,..., p
, where At[i]
is the constant value of a time periodic matrix A(t)
of period T′
for t ∈ [ts[i],ts[i+1])
, if i < p
, or t ∈ [ts[i],T′)
, if i = p
. It is assumed that At[i] := At[mod(i-1,p)+1]
and ts[i] := ts[mod(i-1,p)+1]
for arbitrary i
. All component matrices must have the same dimensions. The component matrices At
, the switching times ts
, the period T
and the number of subperiods k
can be accessed via A.values
, A.ts
, A.period
, and A.nperiod
, respectively.
Base.size
— Methodsize(A::PM)
size(A::PM[, dim])
Return a tuple of two integers containing the dimensions of the continuous-time periodic matrix A
of type PM
, where PM
is one of the types PeriodicFunctionMatrix
, HarmonicArray
, PeriodicTimeSeriesMatrix
, PeriodicSwitchingMatrix
, PeriodicSymbolicMatrix
or PeriodicFunctionMatrix
. Optionally you can specify a dimension dim
to just get the length of that dimension.
Base.eltype
— Methodeltype(A::PM)
Determine the type of the elements of the continuous-time periodic matrix A
of type PM
, where PM
is one of the types PeriodicFunctionMatrix
, HarmonicArray
, PeriodicTimeSeriesMatrix
, PeriodicSwitchingMatrix
, PeriodicSymbolicMatrix
or PeriodicFunctionMatrix
.
Base.getindex
— Methodgetindex(A::PM, ind1, ind2)
Return the continuous-time periodic matrix built from the selected ranges [ind1,ind2]
of the elements of the continuous-time periodic matrix A
of type PM
, where PM
is one of the types PeriodicFunctionMatrix
, HarmonicArray
, PeriodicTimeSeriesMatrix
, PeriodicSwitchingMatrix
, PeriodicSymbolicMatrix
or PeriodicFunctionMatrix
. ind1
and ind2
may be integers, integer ranges or colons.
Base.lastindex
— Methodlastindex(A::PM,dim)
Return the last index along dimension dim
of the continuous-time periodic matrix A
of type PM
, where PM
is one of the types PeriodicFunctionMatrix
, HarmonicArray
, PeriodicTimeSeriesMatrix
, PeriodicSwitchingMatrix
, PeriodicSymbolicMatrix
or PeriodicFunctionMatrix
.