Periodic matrix conversions
convert
Conversion between discrete-time and between continuous-time periodic matrix representations.ts2hr
Conversion of a periodic time series matrix to a harmonic array approximation.pfm2hr
Conversion of a periodic function matrix to a harmonic array representation.ts2pfm
Conversion of an interpolated periodic time series matrix to a periodic function matrix.hr2psm
Conversion of a harmonic array representation to a symbolic matrix.psm2hr
Conversion of a periodic symbolic matrix into a harmonic array representation.pm2pa
Conversion of a discrete-time periodic matrix object to a periodic array object.ffm2hr
Conversion of a Fourier function matrix to a harmonic array representation.ffm2psm
Conversion of a Fourier function matrix to a symbolic matrix.hr2bt
Building a block Toeplitz matrix approximation of a harmonic (Fourier) array representation.hr2btupd
Building an updated block Toeplitz matrix approximation of a harmonic (Fourier) array representation.
Base.convert
— Methodconvert(PM1,A::PM2) -> B::PM1
Convert the discrete-time periodic matrix A
of type PM2
to the discrete-time periodic matrix B
of type PM1
, where PM1
and PM2
are of types PeriodicMatrix
, SwitchingPeriodicMatrix
, PeriodicArray
or SwitchingPeriodicArray
.
Base.convert
— Methodconvert(PM1,A::PM2) -> B::PM1
Convert the continuous-time periodic matrix A
of type PM2
to the continuous-time periodic matrix B
of type PM1
, where PM1
and PM2
are of types PeriodicFunctionMatrix
, HarmonicArray
, PeriodicTimeSeriesMatrix
, PeriodicSwitchingMatrix
, PeriodicSymbolicMatrix
or PeriodicFunctionMatrix
.
PeriodicMatrices.ts2hr
— Function ts2hr(A::PeriodicTimeSeriesMatrix; atol = 0, rtol, n, squeeze = true) -> Ahr::HarmonicArray
Compute the harmonic (Fourier) approximation of a periodic matrix specified by a time series matrix. The periodic matrix A(t)
is specified as a continuous-time periodic time series matrix A
, with m
matrices contained in the vector of matrices A.values
, where A.values[k]
is the value of A(t)
at time moment (k-1)T/m
, with T = A.period
being the period. The resulting harmonic approximation Ahr(t)
of A(t)
has the form
p
Ahr(t) = A_0 + ∑ ( Ac_i*cos(i*t*2*π/T)+As_i*sin(i*2*π*t/T) )
i=1
where A_0
is the constant term (the mean value), Ac_i
and As_i
are the coefficient matrices of the i
-th cosinus and sinus terms, respectively. The order of the approximation p
is determined using the maximum order specified by n
(default: n = (m-1)/2
) and the absolute and relative tolerances atol
and rtol
, as follows: p
is the minimum between n
, (m-1)/2
and the maximum index k
such that Ac_k
and/or As_k
are nonzero. The tolerance used to assess nonzero elements is tol = max(atol, rtol*maxnorm)
, where maxnorm
is the maximum norm of the matrices contained in A.values
. The default values of tolerances are atol = 0
and rtol = 10*p*ϵ
, where ϵ
is the working machine precision.
The resulting harmonic approximation Ahr(t)
is returned in the harmonic array object Ahr
(see HarmonicArray
).
PeriodicMatrices.pfm2hr
— Function pfm2hr(A::PeriodicFunctionMatrix; nsample, NyquistFreq) -> Ahr::HarmonicArray
Convert a periodic function matrix into a harmonic array representation. If A(t)
is a periodic function matrix of period T
, then the harmonic array representation Ahr
is determined by applying the fast Fourier transform to the sampled arrays A(iΔ)
, i = 0, ..., k
, where Δ = T/k
is the sampling period and k
is the number of samples specified by the keyword argument nsample = k
(default: k = 128
). If the Nyquist frequency f
is specified via the keyword argument NyquistFreq = f
, then k
is chosen k = 2*f*T
to avoid signal aliasing.
PeriodicMatrices.ts2pfm
— Function ts2pfm(At::PeriodicTimeSeriesMatrix; method = "linear") -> A::PeriodicFunctionMatrix
Compute the periodic function matrix corresponding to an interpolated periodic time series matrix. For the given periodic time series matrix At
, a periodic function matrix A(t)
is defined as the mapping A(t) = t -> etpf(t)
, where etpf(t)
is a periodic interpolation/extrapolation object, as provided in the Interpolations.jl
package. The keyword parameter method
specifies the interpolation/extrapolation method to be used as follows:
method = "constant"
- use periodic B-splines of degree 0 (periodic constant interpolation);
method = "linear"
- use periodic B-splines of degree 1 (periodic linear interpolation) (default);
method = "quadratic"
- use periodic B-splines of degree 2 (periodic quadratic interpolation);
method = "cubic"
- use periodic B-splines of degree 3 (periodic cubic interpolation).
PeriodicMatrices.hr2psm
— Function hr2psm(Ahr::HarmonicArray, nrange) -> A::Matrix{Num}
Convert a range of harmonic components nrange
of the harmonic array Ahr
to a symbolic matrix A
. The default range is nrange = 0:n
, where n
is the order of the maximum harmonics.
PeriodicMatrices.psm2hr
— Function psm2hr(A::PeriodicSymbolicMatrix; nsample, NyquistFreq) -> Ahr::HarmonicArray
Convert a periodic symbolic matrix into a harmonic array representation. If A(t)
is a periodic symbolic matrix of period T
, then the harmonic array representation Ahr
is determined by applying the fast Fourier transform to the sampled arrays A(iΔ)
, i = 0, ..., k
, where Δ = T/k
is the sampling period and k
is the number of samples specified by the keyword argument nsample = k
(default: k = 128
). If the Nyquist frequency f
is specified via the keyword argument NyquistFreq = f
, then k
is chosen k = 2*f*T
to avoid signal aliasing.
PeriodicMatrices.pm2pa
— Function pm2pa(At::PeriodicMatrix) -> A::PeriodicArray
Convert a discrete-time periodic matrix object into a discrete-time periodic array object.
The discrete-time periodic matrix object At
contains a p
-vector At
of real matrices At[i]
, i = 1,..., p
, the associated time period T
and the number of subperiods k
. The resulting discrete-time periodic array object A
of period T
and number of subperiods k
is built from a m×n×p
real array A
, such that A[:,:,i]
contains At[i]
, for i = 1,..., p
. For non-constant dimensions, the resulting m
and n
are the maximum row and column dimensions, respectively, and the resulting component matrices A[:,:,i]
contain At[i]
, appropriately padded with zero entries.
PeriodicMatrices.ffm2hr
— Function ffm2hr(Afun::FourierFunctionMatrix; atol = 0, rtol = √ϵ, squeeze = true) -> Ahr::HarmonicArray
Compute the harmonic (Fourier) representation of a Fourier periodic matrix object.
The Fourier function matrix object Afun
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 harmonic array object Ahr
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 p′
is the maximum index j
, such that Ac_j
and/or As_j
are nonzero. The tolerance used to assess nonzero elements is tol = max(atol, rtol*maxnorm)
, where maxnorm
is the maximum absolute value of the coefficients ac_i
and as_i
in Afun(t)
. The default values of tolerances are atol = 0
and rtol = √ϵ
, where ϵ
is the working machine precision. The resulting harmonic approximation Ahr(t)
is returned in the harmonic array object Ahr
(see HarmonicArray
).
PeriodicMatrices.ffm2psm
— Function ffm2psm(Af::FourierFunctionMatrix, nrange atol = 0, rtol = 10*n*ϵ,) -> A::Matrix{Num}
Convert a range of harmonic components nrange
of the Fourier function matrix Af
to a symbolic matrix A
. The default range is nrange = 0:n
, where n
is the order of the maximum harmonics. The tolerance used to assess nonzero coefficients is tol = max(atol, rtol*maxnorm)
, where maxnorm
is the maximum absolute value of the coefficients ac_i
and as_i
in Af(t)
. The default values of tolerances are atol = 0
and rtol = 10*n*ϵ
, where ϵ
is the working machine precision.
PeriodicMatrices.hr2bt
— Function hr2bt(Ahr::HarmonicArray, N; P, nperiod]) -> Abt::Matrix
Build the block Toeplitz matrix of a harmonic (Fourier) array representation of a periodic matrix.
The harmonic representation object Ahr
of period T
of a periodic matrix Ahr(t)
of subperiod T′ = T/k
is 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. Ahr(t)
can be equivalently expressed in the Fourier series representation form
p
Ahr(t) = ∑ A_i*exp(im*i*ωh*t) ,
i=-p
where ωh = 2π/T′
, A_i = (Ac_i-im*As_i)/2
and A_{-i} = (Ac_i+im*As_i)/2
. N
is the number of selected harmonic components (or Fourier modes) used for approximation. The keyword parameter P
is the number of full periods to be considered (default: P = 1
) and nperiod
is the number of subperiods to be considered, such that 1 ≤ nperiod ≤ k
(default: nperiod = k
).
For a given number N ≥ p
, if the number of period is P = 1
and the number of subperiods is nperiod = 1
, then the banded block Toeplitz matrix Abt
with (2N+1)×(2N+1)
blocks is constructed
( A_0 A_{-1} … A_{-p} 0 )
( A_1 A_0 ⋱ )
( ⋮ ⋱ ⋱ )
Abt = ( A_p ⋱ A_{-p} )
( ⋱ ⋱ ⋮ )
( 0 A_p … A_0 )
If N < p
, then a truncated full block Toeplitz matrix is built using the first N
harmonic components.
Generally, for given P ≥ 1
and nperiod ≥ 1
, the block Toeplitz matrix Abt
is constructed with (2N*np+1)×(2N*np+1)
blocks, with np = P*nperiod
, such that each A_i
is preceeded in its column by np-1
zero blocks, each A_{-i}
is preceeded in its row by np-1
zero blocks and all diagonal blocks are equal toA_0
.
PeriodicMatrices.hr2btupd
— Function hr2btupd(Ahr::HarmonicArray, N; P, nperiod, shift]) -> Ab::Matrix
Build the updated block Toeplitz matrix of a harmonic (Fourier) array representation of a periodic matrix.
If BT
is the block Toeplitz matrix of the harmonic array representation of the n × n
periodic matrix Ahr
of subperiod T′
(see HarmonicArray
) as constructed with hr2bt
, then the updated matrix Ab = BT-NT is constructed, with NT
a block-diagonal matrix with n × n
diagonal blocks. The k
-th diagonal block of NT
is the diagonal matrix im*(μ + k*ωh)*I
, where μ
is a shift specified via the keyword parameter shift = μ
(default: μ = 0
) and ωh
is the base frequency computed as ωh = 2π*nperiod/(P*T′)
. The value of shift must satisfy 0 ≤ μ ≤ ωh/2
.