API

Main analysis functions

TransitionsInTimeseries.estimate_changesFunction
estimate_changes(config::ChangesConfig, x [,t]) → result

Estimate possible transitions for input timeseries x using the approach specified in the configuration type config, see ChangesConfig for possibilities. t is the time vector corresponding to x, which defaults to eachindex(x).

Return the output as subtype of ChangesResults. The particular form of the output depends on the config and is described in its docstring. Regardless of type, result can always be given to significant_transitions to deduce which possible transitions are statistically significant using a variety of significance tests.

source

Sliding window

TransitionsInTimeseries.SlidingWindowConfigType
SlidingWindowConfig <: ChangesConfig
SlidingWindowConfig(indicators, change_metrics; kwargs...)

A configuration that can be given to estimate_changes. It estimates transitions by a sliding window approach:

  1. Estimate the timeseries of an indicator by sliding a window over the input timeseries.
  2. Estimate changes of an indicator by sliding a window of the change metric over the indicator timeseries.

Both indicators and change metrics are generic Julia functions that input an x::AbstractVector and output an s::Real. Any function may be given and see making custom indicators/change metrics in the documentation for more information on possible optimizations.

indicators can be a single function or a tuple of indicators. Similarly, change_metrics can be a tuple or a single function. If tuples, the length of indicators and change_metrics must match. This way the analysis can be efficiently repeated for many indicators and/or change metrics.

The results output corresponding to SlidingWindowConfig is SlidingWindowResults.

Step 1. is skipped if nothing is provided as indicators, in which case the change metrics are estimated directly from input data.

Keyword arguments

  • width_ind::Int=100, stride_ind::Int=1: width and stride given to WindowViewer to compute the indicator from the input timeseries.
  • width_cha::Int=50, stride_cha::Int=1: width and stride given to WindowViewer to compute the change metric timeseries from the indicator timeseries.
  • whichtime = midpoint: The time vector corresponding to the indicators / change metric timeseries is obtained from t in estimate_changes using the keyword whichtime. Options include:
    • last: use the last timepoint of each window
    • midpoint: use the mid timepoint of each time window
    • first: use first timepoint of each window
    In fact, the indicators time vector is computed simply via
    t_indicator = windowmap(whichtime, t; width_ind, stride_ind)
    t_change = windowmap(whichtime, t_indicator; width_cha, stride_cha)
    so any other function of the time window may be given to extract the time point itself, such as mean or median.
  • T = Float64: Element type of input timeseries to initialize some computations.
source
TransitionsInTimeseries.SlidingWindowResultsType
SlidingWindowResults <: ChangesResults

A struct containing the output of estimate_changes used with SlidingWindowConfig. It can be used for further analysis, visualization, or given to significant_transitions.

It has the following fields that the user may access

  • x: the input timeseries.

  • t: the time vector of the input timeseries.

  • x_indicator, the indicator timeseries (matrix with each column one indicator).

  • t_indicator, the time vector of the indicator timeseries.

  • x_change, the change metric timeseries (matrix with each column one change metric).

  • t_change, the time vector of the change metric timeseries.

  • config::SlidingWindowConfig, what was used for the analysis.

source

Segmented window

TransitionsInTimeseries.SegmentedWindowConfigType
SegmentedWindowConfig <: IndicatorsChangeConfig
SegmentedWindowConfig(indicators, change_metrics, tseg_start, tseg_end; kwargs...)

A configuration that can be given to estimate_changes. It estimates transitions by estimating indicators and changes in user-defined window segments as follows:

  1. For each segment specified, estimate the corresponding indicator timeseries by sliding a window over the input timeseries (within the window segment).
  2. For each segment of the indicator timeseries, estimate a scalar change metric by applying the change metric over the full segment of the indicator timeseries.d

tseg_start, tseg_end are the starts and ends of the window segments (the window segments may overlap, that's okay). indicators, change_metrics are identical as in SlidingWindowConfig.

The results output corresponding to SlidingWindowConfig is SegmentedWindowResults.

Keyword arguments

  • width_ind::Int=100, stride_ind::Int=1, whichtime = midpoint, T = Float64: keywords identical as in SlidingWindowConfig.
  • min_width_cha::Int=typemax(Int): minimal width required to perform the change metric estimation. If a segment is not sufficiently long, the change metric is NaN.
source
TransitionsInTimeseries.SegmentedWindowResultsType
SegmentedWindowResults <: ChangesResults

A struct containing the output of estimate_changes used with SegmentedWindowConfig. It can be used for further analysis, visualization, or given to significant_transitions.

It has the following fields that the user may access

  • x: the input timeseries.

  • t: the time vector of the input timeseries.

  • x_indicator::Vector{Matrix}, with x_indicator[k] the indicator timeseries (matrix with each column one indicator) of the k-th segment.

  • t_indicator::Vector{Vector}, with t_indicator[k] the time vector of the indicator timeseries for the k-th segment.

  • x_change::Matrix, the change metric values with x[k, i] the change metric of the i-th indicator for the k-th segment.

  • t_change, the time vector of the change metric.

  • config::SegmentedWindowConfig, what was used for the analysis.

  • i1::Vector{Int} indices corresponding to start time of each segment.

  • i2::Vector{Int} indices corresponding to end time of each segment.

  • precomp_change_metrics vector containing the precomputed change metrics of each segment.

source

Slope change

TransitionsInTimeseries.SlopeChangeConfigType
SlopeChangeConfig <: ChangesConfig
SlopeChangeConfig(; indicator = nothing, kw...)

A configuration that can be given to estimate_changes. It estimates a change of slope in the timeseries by fitting two connected linear segments to the timeseries, returning the results (i.e., the two-linear fits) as SlopeChangeResults.

Keyword arguments

  • indicator = nothing: if not nothing. Otherwise it should be a function f(x) -> Real. The slope fitting is then done over an indicator of the timeseries, which itself is estimated via a sliding window exactly as in SlidingWindowConfig.
  • width_ind, stride_ind, whichtime: exactly as in SlidingWindowConfig if indicator is not nothing.
source
TransitionsInTimeseries.SlopeChangeResultsType
SlopeChangeResults <: ChangesResults

A struct containing the output of estimate_changes used with SlopeChangeConfig. It can be used for further analysis, visualization, or given to significant_transitions. The only significance type that you can use this with significant_transitions is SlopeChangeSignificance.

It has the following fields that the user may access:

  • x: the input timeseries.
  • t: the time vector of the input timeseries.
  • x_indicator, the indicator timeseries.
  • t_indicator, the time vector of the indicator timeseries.
  • t_change, the time the slope changes.
  • fitparams = a, b, c, d, the fitted linear coefficients, a + b*t before. t_change and c + d*t after t_change.
source

Significance testing

TransitionsInTimeseries.significant_transitionsFunction
significant_transitions(res::ChangesResults, signif::Significance)

Estimate significant transtions in res using the method described by signif. Return flags, a Boolean matrix with identical size as the changes stored in res (which typically is stored in the field res.x_change). flags is true wherever a change metric of res is deemed significant.

source
TransitionsInTimeseries.SurrogatesSignificanceType
SurrogatesSignificance <: Significance
SurrogatesSignificance(; kwargs...)

A configuration struct for significance testing significant_transitions using timeseries surrogates.

Keyword arguments

  • surromethod = RandomFourier(): method to generate surrogates
  • n = 1000: how many surrogates to generate
  • rng = Random.default_rng(): random number generator for the surrogates
  • p = 0.05: threshold for significance of the p-value
  • tail = :both: tail type used, see below

Description

When used with ChangesResults, significance is estimated as follows: n surrogates from the input timeseries are generated using surromethod, which is any Surrogate subtype provided by TimeseriesSurrogates.jl. For each surrogate, the indicator and then change metric timeseries is extracted. The values of the surrogate change metrics form a distribution of values (one at each time point). The value of the original change metric is compared to that of the surrogate distribution and a p-value is extracted according to the specified tail. The p-value is compared with p to claim significance. After using SurrogatesSignificance, you may access the full p-values before thresholding in the field .pvalues (to e.g., threshold with different p).

The p-value is simply the proportion of surrogate change metric values that exceed (for tail = :right) or subseed (tail = :left) the original change metric at each given time point. Use tail = :left if the surrogate data are expected to have higher change metric values. This is the case for statistics that quantify entropy. For statistics that quantify autocorrelation, use tail = :right instead. For anything else, use the default tail = :both. An iterable of tail values can also be given, in which case a specific tail is used for each change metric in ChangesResults.

Note that the raw p-values can be accessed in the field .pvalues, after calling the significant_transitions function with SurrogatesSignificance, in case you wish to obtain a different threshold of the p-values.

source
TransitionsInTimeseries.ThresholdSignificanceType
ThresholdSignificance(threshold::Real; tail = :right) <: Significance

A configuration struct for significance testing in significant_transitions. Significance is estimated by comparing the value of each change metric with the given threshold. Values that exceed the threshold (if tail = :right) or subseed the threshold (if tail = :left) are deemed significant. If tail = :both then either condition is checked.

source
TransitionsInTimeseries.SigmaSignificanceType
SigmaSignificance(; factor = 3.0, tail = :both) <: Significance

A configuration struct for significance testing significant_transitions. When used with ChangesResults, significance is estimated by comparing how many standard deviations (σ) the value exceeds the mean value (μ). Values that exceed (if tail = :right) μ + factor*σ, or subseed (if tail = :left) μ - factor*σ are deemed significant. If tail = :both then either condition is checked.

factor can also be a vector of values, in which case a different value is used for each change metric.

See also QuantileSignificance.

source
TransitionsInTimeseries.QuantileSignificanceType
QuantileSignificance(; p = 0.95, tail = :both) <: Significance

A configuration struct for significance testing significant_transitions. When used with ChangesResults, significance is estimated by comparing the value of each change metric with its p-quantile. Values that exceed the p-quantile (if tail = :right) or subseed the 1-p-quantile (if tail = :left) are deemed significant. If tail = :both then either condition is checked.

QuantileSignficance guarantees that some values will be significant by the very definition of what a quantile is. See also SigmaSignificance that is similar but does not have this guarantee.

source
TransitionsInTimeseries.SlopeChangeSignificanceType
SlopeChangeSignificance(; moe_slope, moe_offset, slope_diff = moe_slope, pvalue = 0.05)

Test whether the result of SlopeChangeResults is statistically significant.

Two tests are done:

  1. Check whether the margin of error of the fitted parameters a, b, c, d of the two linear segments a + b*t, c + d*t is less than the specified margins of error, for a chosen pvalue.
  2. Test that the two slopes b, d have difference greater than slope_diff.

The Boolean & of the above two is the final test.

The margin of error is simply half the size of the confidence interval, also known as radius of the confidence interval.

source

Indicators

Note that any Julia function can be an indicator or change metric, so the list here is only just a couple of indicators directly implemented in this package.

Value distribution

Statistics.meanFunction
mean(A::AbstractArray, w::AbstractWeights[, dims::Int])

Compute the weighted mean of array A with weight vector w (of type AbstractWeights). If dim is provided, compute the weighted mean along dimension dims.

Examples

n = 20
x = rand(n)
w = rand(n)
mean(x, weights(w))
source
StatsBase.skewnessFunction
skewness(v, [wv::AbstractWeights], m=mean(v))

Compute the standardized skewness of a real-valued array v, optionally specifying a weighting vector wv and a center m.

source
StatsBase.kurtosisFunction
kurtosis(v, [wv::AbstractWeights], m=mean(v))

Compute the excess kurtosis of a real-valued array v, optionally specifying a weighting vector wv and a center m.

source

Critical Slowing Down

Statistics.varFunction
var(x::AbstractArray, w::AbstractWeights, [dim]; mean=nothing, corrected=false)

Compute the variance of a real-valued array x, optionally over a dimension dim. Observations in x are weighted using weight vector w. The uncorrected (when corrected=false) sample variance is defined as:

\[\frac{1}{\sum{w}} \sum_{i=1}^n {w_i\left({x_i - μ}\right)^2 }\]

where $n$ is the length of the input and $μ$ is the mean. The unbiased estimate (when corrected=true) of the population variance is computed by replacing $\frac{1}{\sum{w}}$ with a factor dependent on the type of weights used:

  • AnalyticWeights: $\frac{1}{\sum w - \sum {w^2} / \sum w}$
  • FrequencyWeights: $\frac{1}{\sum{w} - 1}$
  • ProbabilityWeights: $\frac{n}{(n - 1) \sum w}$ where $n$ equals count(!iszero, w)
  • Weights: ArgumentError (bias correction not supported)
source
var(ce::CovarianceEstimator, x::AbstractVector; mean=nothing)

Compute the variance of the vector x using the estimator ce.

source
TransitionsInTimeseries.ar1_whitenoiseFunction
ar1_whitenoise(x::AbstractVector)

Return the AR1 regression coefficient of a time series x by computing the analytic solution of the least-square parameter estimation under white-noise assumption for the data-generating process:

\[\theta = \sum_{i=2}^{n} x_i x_{i-1} / \sum_{i=2}^{n} x_{i-1}^2\]

source

Spectrum

TransitionsInTimeseries.LowfreqPowerSpectrumType
LowfreqPowerSpectrum(; q_lofreq = 0.1)

Return a PrecomputableFunction containing all the necessary fields to generate a PrecomputedLowfreqPowerSpectrum. The latter can be initialized by precompute:

lfps = precompute( LowfreqPowerSpectrum() )

Keyword arguments:

  • q_lofreq: a number between 0 and 1 that characterises which portion of the

frequency spectrum is considered to be low. For instance, q_lofreq = 0.1 implies that the lowest 10% of frequencies are considered to be the low ones.

source
TransitionsInTimeseries.PrecomputedLowfreqPowerSpectrumType
PrecomputedLowfreqPowerSpectrum

A struct containing all the precomputed fields to efficiently perform repetitive computation of the low-frequency power spectrum (LFPS), a number between 0 and 1 that characterizes the amount of power contained in the low frequencies of the power density spectrum of x. Once lfps::PrecomputedLowfreqPowerSpectrum is initialized, it can be used as a function to obtain the LFPS of x::AbstractVector by:

lfps(x)
source

Nonlinear dynamics

Indicators that come from nonlinear timeseries analysis typically quantify some entropy-like or complexity measure from the timeseries. Thousands (literally) such measures are provided out of the box by the ComplexityMeasures.jl package. Given that any of these may be used as an indicator or change metric, we made the decision to not copy-paste any measure here, as it is easy for the user to use any of them.

For example, using the permutation entropy as an indicator is as simple as doing

using ComplexityMeasures
est = OrdinalPatterns(; m = 3) # order 3
# create a function that given timeseries returns permutation entropy
indicator = x -> entropy_normalized(est, x)

and giving the created indicator to e.g., SlidingWindowConfig.

Change metrics

Slope

TransitionsInTimeseries.kendalltauFunction
kendalltau(x::AbstractVector)

Compute the kendall-τ correlation coefficient of the time series x. kendalltau can be used as a change metric focused on trend.

source
TransitionsInTimeseries.spearmanFunction
spearman(x::AbstractVector)

Compute the spearman correlation coefficient of the time series x. spearman can be used as a change metric focused on trend.

source

Value distribution differences

TransitionsInTimeseries.difference_of_meansFunction
difference_of_means(x::AbstractArray)

Return the absolute difference of the means of the first and second halfs of x. difference_of_means can be used as a change metric focused on value differences. Creating similar statistical differences using other moments instead of mean is trivial. In fact, the source of difference_of_means is just:

# assumes 1-based indexing
n = length(x)
x1 = view(x, 1:n÷2)
x2 = view(x, (n÷2 + 1):n)
return abs(mean(x1) - mean(x2))

difference_of_means can also sensibly be used for windows of size 2, in which case the change metric timeseries is the same as the abs.(diff(...)) of the indicator timeseries.

source

Make your own indicator/metric!

The only difference between what is an "indicator" and what is a "change metric" is purely conceptual. As far as the code base of TransitionsInTimeseries.jl is concerned, they are both functions f: x::AbstractVector{Real} -> f(x)::Real. As a user you may give any such function for an indicator or change metric.

There are situations where you may optimize such a function based on knowledge of input x type and length, in which case you may use PrecomputableFunction:

Surrogates

For the surrogate generation, you can use any subtype of Surrogate defined in Timeseriessurrogates.jl.

Sliding windows

TransitionsInTimeseries.WindowViewerType
WindowViewer(x; width, stride)

Initialize an iterator that generates views over the given timeseries x based on a window with a given width, incrementing the window views with the given stride. You can use this directly with map, such as map(std, WindowViewer(x, ...)) would give you the moving-window-timeseries of the std of x.

If not given, the keywords width, stride are respectively taken as default_window_width(x) and 1.

source
TransitionsInTimeseries.windowmapFunction
windowmap(f::Function, x::AbstractVector; kwargs...) → mapped_f

A shortcut for first generating a wv = WindowViewer(x; kwargs...) and then applying mapped_f = map(f, wv). If x is accompanied by a time vector t, you probably also want to call this function with t instead of x and with one of mean, midpoint, midvalue as f to obtain a time vector for the mapped_f output.

source

Load data

Visualization

TransitionsInTimeseries.plot_indicator_changesFunction
plot_indicator_changes(res) → (fig, axs)

Return fig::Figure and axs::Matrix{Axis}, on which res::ChangesResults has been visualised.

Keyword arguments:

  • colors = default_colors sets the colors of the line plots that are to

be cycled through. The default correspond to the color scheme of Julia Dynamics.

  • `indicatornames = defaultindicatorlabel(res), chametricnames =

defaultchametriclabel(res)sets the labels for the indicators and the change metrics, with the default inferring them from the names ofres.config.indicatorsandres.config.change_metrics`.

  • accent_linewidth = 3 sets the line width for the original signals (the

surrogates have linewidth = 1)

source
TransitionsInTimeseries.plot_significance!Function
plot_significance!(axs, res, signif)

Update the axs::Matrix{Axis} of a figure obtained with plot_indicator_changes(res) with the signif::SurrogatesSignificance.

Keyword arguments:

  • flags = nothing provides the significance flags, for instance obtained by

thresholding the pvalues.

  • nsurro = 20 sets the number of surrogates to plot.
source

Utils

default_window_width