Annex functions
Functions
This package also includes the following general-purpose functions useful for time series analysis.
Building blocks for time series algorithms
Base
MessyTimeSeries.pushcopy!
— Functionpushcopy!(collection::AbstractVector, item::AbstractArray)
Push copy of item
into collection
.
Convenient mathematical and statistical operations
MessyTimeSeries.no_combinations
— Functionno_combinations(n::Int64, k::Int64)
Compute the binomial coefficient of n
observations and k
groups, for big integers.
Examples
julia> no_combinations(1000000,100000)
7.333191945934207610471288280331309569215030711272858517142085449641265002716664e+141178
MessyTimeSeries.rand_without_replacement
— Functionrand_without_replacement(rng::StableRNGs.LehmerRNG, nT::Int64, d::Int64)
Draw length(P)-d
elements from the positional vector P
without replacement.
P
is permanently changed in the process.
Examples
julia> rand_without_replacement(StableRNG(1), 20, 5)
5-element Array{Int64,1}:
1
9
14
19
20
rand_without_replacement(rng::StableRNGs.LehmerRNG, n::Int64, T::Int64, d::Int64)
Draw length(P)-d
elements from the positional vector P
without replacement.
In the sampling process, no more than n-1 elements are removed for each point in time. P
is permanently changed in the process.
MessyTimeSeries.soft_thresholding
— Functionsoft_thresholding(z::Float64, ζ::Float64)
Soft thresholding operator.
MessyTimeSeries.solve_discrete_lyapunov
— Functionsolve_discrete_lyapunov(A::AbstractArray{Float64,2}, Q::SymMatrix)
Use a bilinear transformation to convert the discrete Lyapunov equation to a continuous Lyapunov equation, which is then solved using BLAS.
The notation used for representing the discrete Lyapunov equation is
$P - APA' = Q$,
where $P$ and $Q$ are symmetric. This equation is transformed into
$B'P + PB = -C$
References
Kailath (1980, page 180)
Convergence check
MessyTimeSeries.check_bounds
— Functioncheck_bounds(X::Real, LB::Real, UB::Real)
Check whether X
is larger or equal than LB
and lower or equal than UB
check_bounds(X::Real, LB::Real)
Check whether X
is larger or equal than LB
MessyTimeSeries.isconverged
— Functionisconverged(new::Float64, old::Float64, tol::Float64, ε::Float64, increasing::Bool)
Check whether new
is close enough to old
.
Arguments
new
: new objective or lossold
: old objective or losstol
: toleranceε
: small Float64increasing
: true ifnew
increases, at each iteration, with respect toold
Parameter transformations
MessyTimeSeries.get_bounded_log
— Functionget_bounded_log(Θ_unbound::Float64, MIN::Float64)
Compute parameters with bounded support using a generalised log transformation.
MessyTimeSeries.get_unbounded_log
— Functionget_unbounded_log(Θ_bound::Float64, MIN::Float64)
Compute parameters with unbounded support using a generalised log transformation.
MessyTimeSeries.get_bounded_logit
— Functionget_bounded_logit(Θ_unbound::Float64, MIN::Float64, MAX::Float64)
Compute parameters with bounded support using a generalised logit transformation.
MessyTimeSeries.get_unbounded_logit
— Functionget_unbounded_logit(Θ_bound::Float64, MIN::Float64, MAX::Float64)
Compute parameters with unbounded support using a generalised logit transformation.
Sample statistics for incomplete data
MessyTimeSeries.mean_skipmissing
— Functionmean_skipmissing(X::AbstractVector{Float64})
mean_skipmissing(X::AbstractVector{Union{Missing, Float64}})
Compute the mean of the observed values in X
.
Examples
julia> mean_skipmissing([1.0; missing; 3.0])
2.0
mean_skipmissing(X::AbstractMatrix{Float64})
mean_skipmissing(X::AbstractMatrix{Union{Missing, Float64}})
Compute the mean of the observed values in X
column wise.
Examples
julia> mean_skipmissing([1.0 2.0; missing 3.0; 3.0 5.0])
3-element Array{Float64,1}:
1.5
3.0
4.0
MessyTimeSeries.std_skipmissing
— Functionstd_skipmissing(X::AbstractVector{Float64})
std_skipmissing(X::AbstractVector{Union{Missing, Float64}})
Compute the standard deviation of the observed values in X
.
Examples
julia> std_skipmissing([1.0; missing; 3.0])
1.4142135623730951
std_skipmissing(X::AbstractMatrix{Float64})
std_skipmissing(X::AbstractMatrix{Union{Missing, Float64}})
Compute the standard deviation of the observed values in X
column wise.
Examples
julia> std_skipmissing([1.0 2.0; missing 3.0; 3.0 5.0])
3-element Array{Float64,1}:
0.7071067811865476
NaN
1.4142135623730951
MessyTimeSeries.sum_skipmissing
— Functionsum_skipmissing(X::AbstractVector{Float64})
sum_skipmissing(X::AbstractVector{Union{Missing, Float64}})
Compute the sum of the observed values in X
.
Examples
julia> sum_skipmissing([1.0; missing; 3.0])
4.0
sum_skipmissing(X::AbstractMatrix{Float64})
sum_skipmissing(X::AbstractMatrix{Union{Missing, Float64}})
Compute the sum of the observed values in X
column wise.
Examples
julia> sum_skipmissing([1.0 2.0; missing 3.0; 3.0 5.0])
3-element Array{Float64,1}:
3.0
3.0
8.0
Time-series operations
Foundations
MessyTimeSeries.companion_form
— Functioncompanion_form(Ψ::AbstractArray{Float64,2}; extended::Bool=false)
Construct the companion form matrix from the generic coefficients Ψ.
If extended is true, it increases the typical dimension of the companion matrix by n rows.
MessyTimeSeries.lag
— Functionlag(X::FloatArray, p::Int64)
Construct the data required to run a standard vector autoregression.
Arguments
X
: observed measurements (nxT
), wheren
andT
are the number of series and observations.p
: number of lags in the vector autoregression
Output
X_{t}
X_{t-1}
MessyTimeSeries.diff2
— Functiondiff2(A::AbstractArray, dims::Integer)
Return double-differenced data.
MessyTimeSeries.diff_or_diff2
— Functiondiff_or_diff2(A::AbstractArray, dims::Integer, use_diff::Bool)
Use either diff or diff2 depending on the value taken by use_diff
.
Interpolation and moving averages
MessyTimeSeries.centred_moving_average
— Functioncentred_moving_average(X::Union{FloatMatrix, JMatrix{Float64}}, n::Int64, T::Int64, window::Int64)
Compute the centred moving average of X
.
Arguments
X
: observed measurements (nxT
)n
andT
are the number of series and observationswindow
is the total number of observations (lagging, current and leading) included in the average
centred_moving_average(X::Union{FloatMatrix, JMatrix{Float64}}, window::Int64)
Compute the centred moving average of X
.
Arguments
X
: observed measurements (nxT
)window
is the total number of observations (lagging, current and leading) included in the average
MessyTimeSeries.forward_backwards_rw_interpolation
— Functionforward_backwards_rw_interpolation(X::JMatrix{Float64}, n::Int64, T::Int64)
Interpolate each non-stationary series in X
, in turn, using a random walk logic both forward and backwards in time.
Arguments
X
: observed measurements (nxT
)n
andT
are the number of series and observations
forward_backwards_rw_interpolation(X::FloatMatrix, n::Int64, T::Int64)
Return X
.
Arguments
X
: observed measurements (nxT
)n
andT
are the number of series and observations
forward_backwards_rw_interpolation(X::JMatrix{Float64})
Interpolate each non-stationary series in X
, in turn, using a random walk logic both forward and backwards in time.
Arguments
X
: observed measurements (nxT
)
MessyTimeSeries.interpolate_series
— Functioninterpolate_series(X::JMatrix{Float64}, n::Int64, T::Int64)
Interpolate each series in X
, in turn, by replacing missing observations with the sample average of the non-missing values.
Arguments
X
: observed measurements (nxT
)n
andT
are the number of series and observations
interpolate_series(X::FloatMatrix, n::Int64, T::Int64)
Return X
.
Arguments
X
: observed measurements (nxT
)n
andT
are the number of series and observations
interpolate_series(X::Union{FloatMatrix, JMatrix{Float64}})
Interpolate each series in X
, in turn, by replacing missing observations with the sample average of the non-missing values.
Arguments
X
: observed measurements (nxT
)
Standardisation
MessyTimeSeries.demean
— Functiondemean(X::FloatVector)
demean(X::FloatMatrix)
demean(X::JVector{Float64})
demean(X::JMatrix{Float64})
Demean X
.
Examples
julia> demean([1.0; 1.5; 2.0; 2.5; 3.0])
5-element Array{Float64,1}:
-1.0
-0.5
0.0
0.5
1.0
julia> demean([1.0 3.5 1.5 4.0 2.0; 4.5 2.5 5.0 3.0 5.5])
2×5 Array{Float64,2}:
-1.4 1.1 -0.9 1.6 -0.4
0.4 -1.6 0.9 -1.1 1.4
MessyTimeSeries.standardise
— Functionstandardise(X::FloatVector)
standardise(X::FloatMatrix)
standardise(X::JVector{Float64})
standardise(X::JMatrix{Float64})
Standardise X
.
Examples
julia> standardise([1.0; 1.5; 2.0; 2.5; 3.0])
5-element Array{Float64,1}:
-1.2649110640673518
-0.6324555320336759
0.0
0.6324555320336759
1.2649110640673518
julia> standardise([1.0 3.5 1.5 4.0 2.0; 4.5 2.5 5.0 3.0 5.5])
2×5 Array{Float64,2}:
-1.08173 0.849934 -0.695401 1.23627 -0.309067
0.309067 -1.23627 0.695401 -0.849934 1.08173
Index
MessyTimeSeries.centred_moving_average
MessyTimeSeries.check_bounds
MessyTimeSeries.companion_form
MessyTimeSeries.demean
MessyTimeSeries.diff2
MessyTimeSeries.diff_or_diff2
MessyTimeSeries.forward_backwards_rw_interpolation
MessyTimeSeries.get_bounded_log
MessyTimeSeries.get_bounded_logit
MessyTimeSeries.get_unbounded_log
MessyTimeSeries.get_unbounded_logit
MessyTimeSeries.interpolate_series
MessyTimeSeries.isconverged
MessyTimeSeries.lag
MessyTimeSeries.mean_skipmissing
MessyTimeSeries.no_combinations
MessyTimeSeries.pushcopy!
MessyTimeSeries.rand_without_replacement
MessyTimeSeries.soft_thresholding
MessyTimeSeries.solve_discrete_lyapunov
MessyTimeSeries.standardise
MessyTimeSeries.std_skipmissing
MessyTimeSeries.sum_skipmissing