timeseries.recurrence_plot

Provides classes for the analysis of dynamical systems and time series based on recurrence plots, including measures of recurrence quantification analysis (RQA) and recurrence network analysis.

class pyunicorn.timeseries.recurrence_plot.RecurrencePlot(time_series: ndarray[Any, dtype[_ScalarType_co]], metric: str = 'supremum', normalize: bool = False, missing_values: bool = False, sparse_rqa: bool = False, silence_level: int = 0, **kwargs)[source]

Bases: Cached

Class RecurrencePlot for generating and quantitatively analyzing recurrence plots.

The RecurrencePlot class supports the construction of recurrence plots from multi-dimensional time series, optionally using embedding. Currently, manhattan, euclidean and supremum norms are provided for measuring distances in phase space.

Methods for calculating commonly used measures of recurrence quantification analysis (RQA) are provided, e.g., determinism, maximum diagonal line length and laminarity. The definitions of these measures together with a review of the theory and applications of recurrence plots are given in [Marwan2007].

Examples:

  • Create an instance of RecurrencePlot with a fixed recurrence threshold and without embedding:

    RecurrencePlot(time_series, threshold=0.1)
    
  • Create an instance of RecurrencePlot with a fixed recurrence threshold in units of STD and without embedding:

    RecurrencePlot(time_series, threshold_std=0.03)
    
  • Create an instance of RecurrencePlot at a fixed (global) recurrence rate and using time delay embedding:

    RecurrencePlot(time_series, dim=3, tau=2,
                   recurrence_rate=0.05).recurrence_rate()
    
N: int

The number of state vectors (number of lines and rows) of the RP.

R

The recurrence matrix.

__cache_state__() Tuple[Hashable, ...][source]

Hashable tuple of mutable object attributes, which will determine the instance identity for ALL cached method lookups in this class, in addition to the built-in object id(). Returning an empty tuple amounts to declaring the object immutable in general. Mutable dependencies that are specific to a method should instead be declared via @Cached.method(attrs=(…)).

NOTE: A subclass is responsible for the consistency and cost of this state descriptor. For example, hashing a large array attribute may be circumvented by declaring it as a property, with a custom setter method that increments a dedicated mutation counter.

__init__(time_series: ndarray[Any, dtype[_ScalarType_co]], metric: str = 'supremum', normalize: bool = False, missing_values: bool = False, sparse_rqa: bool = False, silence_level: int = 0, **kwargs)[source]

Initialize an instance of RecurrencePlot.

Either recurrence threshold threshold/threshold_std, recurrence rate recurrence_rate or local recurrence rate local_recurrence_rate have to be given as keyword arguments.

Embedding is only supported for scalar time series. If embedding dimension dim and delay tau are both given as keyword arguments, embedding is applied. Multidimensional time series are processed as is by default.

Attention:

The sparse_rqa feature is experimental and currently only works for fixed threshold and the supremum metric.

Parameters:
  • time_series (2D array (time, dimension)) – The time series to be analyzed, can be scalar or multi-dimensional.

  • metric (str) – The metric for measuring distances in phase space (“manhattan”, “euclidean”, “supremum”).

  • normalize (bool) – Decide whether to normalize the time series to zero mean and unit standard deviation.

  • missing_values (bool) – Toggle special treatment of missing values in RecurrencePlot.time_series.

  • sparse_rqa (bool) – Toggles sequential RQA computation using less memory for use with long time series.

  • skip_recurrence (bool) – Skip calculation of recurrence matrix within RP class (e.g. when overloading respective methods in child class)

  • silence_level (int) – Inverse level of verbosity of the object.

  • threshold (number) – The recurrence threshold keyword for generating the recurrence plot using a fixed threshold.

  • threshold_std (number) – The recurrence threshold keyword for generating the recurrence plot using a fixed threshold in units of the time series’ STD.

  • recurrence_rate (number) – The recurrence rate keyword for generating the recurrence plot using a fixed recurrence rate.

  • local_recurrence_rate (number) – The local recurrence rate keyword for generating the recurrence plot using a fixed local recurrence rate (same number of recurrences for each state vector).

  • adaptive_neighborhood_size (number) – The adaptive neighborhood size parameter for generating recurrence plots based on the algorithm in [Xu2008].

  • dim (number) – The embedding dimension.

  • tau (number) – The embedding delay.

__str__()[source]

Returns a string representation.

average_diaglength(l_min=2, resampled_dist=None)[source]

Return diagonal line-based RQA measure average diagonal line length \(L\).

\(L\) is defined as the average length of diagonal lines (of at least length \(l_min\)).

Parameters:
  • l_min (number) – The minimum diagonal line length.

  • resampled_dist (1D array (integer)) – resampled frequency distribution of diagonal lines

Return number:

the average diagonal line length \(L\).

average_vertlength(v_min=2, resampled_dist=None)[source]

Return vertical line-based RQA measure average vertical line length \(TT\).

\(TT\) is defined as the average vertical line length (of at least length \(v_min\)) and is also called trapping time \(TT\).

Parameters:
  • v_min (number) – The minimal vertical line length.

  • resampled_dist (1D array (integer)) – resampled frequency distribution of vertical lines

Return number:

the trapping time \(TT\).

average_white_vertlength(w_min=1)[source]

Return white vertical line-based RQA measure average white vertical line length.

It is defined as the average white vertical line length (of at least length \(w_min\)) and is also called mean recurrence time.

Reference: [Ngamga2007].

Parameters:

w_min (number) – The minimal white vertical line length.

Return number:

the mean recurrence time.

static bootstrap_distance_matrix(embedding, metric, M)[source]

Return bootstrap samples from distance matrix.

Parameters:
  • embedding (2D array (time, embedding dimension)) – The phase space trajectory.

  • metric (str) – The metric for measuring distances in phase space (“manhattan”, “euclidean”, “supremum”).

  • M (int) – Number of bootstrap samples

Return type:

1D array (“float32”)

Returns:

the bootstrap samples from distance matrix.

complexity_entropy()[source]

Returns the complexity entropy for each dimension of the time series.

Reference: [Ribeiro2011]

Return type:

double

Returns:

Complexity entropy of the embedded time series

determinism(l_min=2, resampled_dist=None)[source]

Return diagonal line-based RQA measure determinism \(DET\).

\(DET\) is defined as the ratio of recurrence points that form diagonal structures (of at least length \(l_min\)) to all recurrence points.

Parameters:
  • l_min (number) – The minimum diagonal line length.

  • resampled_dist (1D array (integer)) – resampled frequency distribution of diagonal lines

Return number:

the determinism \(DET\).

diag_entropy(l_min=2, resampled_dist=None)[source]

Return diagonal line-based RQA measure diagonal line entropy \(ENTR\).

\(ENTR\) is defined as the entropy of the probability to find a diagonal line of exactly length l in the RP - reflects the complexity of the RP with respect to diagonal lines.

Parameters:
  • l_min (number) – The minimal diagonal line length.

  • resampled_dist (1D array (integer)) – resampled frequency distribution of diagonal lines

Return number:

the diagonal line-based entropy \(ENTR\).

diagline_dist()[source]

Return the frequency distribution of diagonal line lengths \(P(l-1)\).

Note that entry \(P(l-1)\) contains the number of diagonal lines of length \(l\). Thus, \(P(0)\) counts lines of length \(1\), \(P(1)\) counts lines of length \(2\), asf. The main diagonal is not counted, hence \(P(N)\) will always be \(0\).

Note

Experimental handling of missing values. Diagonal lines touching lines and blocks of missing entries in the recurrence matrix are not counted.

Return type:

1D array (int32)

Returns:

the frequency distribution of diagonal line lengths \(P(l-1)\).

distance_matrix(metric: str)[source]

Return phase space distance matrix \(D\) according to the chosen metric.

Parameters:

metric (str) – The metric for measuring distances in phase space (“manhattan”, “euclidean”, “supremum”).

Return type:

2D square array

Returns:

the phase space distance matrix \(D\)

static embed_time_series(time_series, dim, tau)[source]

Return a time series’ delay embedding.

Returns a Numpy array containing a delay embedding of the time series using embedding dimension dim and time delay tau.

Parameters:
  • time_series (1D array) – The scalar time series to be embedded.

  • dim (int) – The embedding dimension.

  • tau (int) – The embedding delay.

Return type:

2D array (time, dimension)

Returns:

the embedded phase space trajectory.

property embedding: ndarray

The embedded time series / phase space trajectory (time, embedding dimension).

euclidean_distance_matrix()[source]

Return the euclidean distance matrix from an embedding of a time series.

Return type:

2D square array

Returns:

the euclidean distance matrix.

laminarity(v_min=2, resampled_dist=None)[source]

Return vertical line-based RQA measure laminarity \(LAM\).

\(LAM\) is defined as the ratio of recurrence points that form vertical structures (of at least length \(v_min\)) to all recurrence points.

Parameters:
  • v_min (number) – The minimal vertical line length.

  • resampled_dist (1D array (integer)) – resampled frequency distribution of vertical lines

Return number:

the laminarity \(LAM\).

static legendre_coordinates(x, dim=3, t=None, p=None, tau_w='est')[source]

Return a phase space trajectory reconstructed using orthogonal polynomial filters.

The reconstructed state vector components are the zero-th to (dim-1)-th derivatives of the (possibly irregularly spaced) time series x as estimated by folding with the orthogonal polynomial filters that correspond to the sequence of measurement time points t.

This is a generalization for irregularly spaced time series of the “Legendre coordinates” introduced in Gibson et al. (1992).

Parameters:
  • x (array-like) – Time series values

  • dim (int) – Dimension > 0 of reconstructed phase space. Default: 3

  • t (array-like or None) – Optional array of measurement time points corresponding to the values in x. Default: [0,…,x.size-1]

  • p (int > 0 or None) – No. of past and future time points to use for the estimation. Default: dim or determined by tau_w if given

  • tau_w (float > 0 or "est" or None) – Optional (average) window width to use in determining p when p = None. Following Gibson et al. (1992), this should be about sqrt(3/< x**2 >) * std(x), or about a quarter period. If “est”, this is estimated iteratively, starting with 4 * (max(t)-min(t)) / (N-1) and estimating x’ from that.

Return type:

2D array [observation index, dimension index]

Returns:

Estimated derivatives. Rows are reconstructed state vectors.

manhattan_distance_matrix()[source]

Return the manhattan distance matrix from an embedding of a time series.

Return type:

2D square array

Returns:

the manhattan distance matrix.

max_diaglength()[source]

Return diagonal line-based RQA measure maximum diagonal line length \(L_max\).

\(L_max\) is defined as the maximal length of a diagonal line in the recurrence matrix.

Return number:

the maximal diagonal line length \(L_max\).

max_vertlength()[source]

Return vertical line-based RQA measure maximal vertical line length \(V_max\).

\(V_max\) is defined as the maximal length of a vertical line of the recurrence matrix.

Return number:

the maximal vertical line length \(V_max\).

max_white_vertlength()[source]

Return white vertical line-based RQA measure maximal white vertical line length.

It is defined as the maximal length of a white vertical line of the recurrence matrix and corresponds to the maximum recurrence time occuring in the time series.

Return number:

the maximal white vertical line length.

mean_recurrence_time(w_min=1)[source]

Alias for average_white_vertlength() (see description there).

metric

The metric used for measuring distances in phase space.

missing_values

Controls special treatment of missing values in RecurrencePlot.time_series.

static normalize_time_series(time_series)[source]

Normalize each component of a time series in place.

Works also for complex valued time series.

Note

Modifies the given array in place!

Parameters:

time_series (2D array (time, dimension)) – The time series to be normalized.

permutation_entropy(normalize=True)[source]

Returns the permutation entropy for an embedded time series. An embedding of 3 <= embedding dimension <= 7 is recommended for this method.

Reference: [Bandt2002]

Return type:

double

Returns:

Permutation entropy of the embedded time series

recurrence_matrix()[source]

Return the current recurrence matrix \(R\).

Return type:

2D square Numpy array

Returns:

the current recurrence matrix \(R\).

recurrence_probability(lag=0)[source]

Return the recurrence probability. This is the probability, that the trajectory is recurrent after ‘lag’ time steps.

Contributed by Jan H. Feldhoff.

Return number:

the recurrence probability

recurrence_rate()[source]

Return the recurrence rate \(RR\).

RR gives the percentage of black dots in the recurrence plot.

Return number:

the recurrence rate \(RR\).

static rejection_sampling(dist, M)[source]

Rejection sampling of discrete frequency distribution.

Use simple rejection sampling algorithm for computing a resampled version of a given frequency distribution with discrete support.

Parameters:
  • dist (1D array (integer)) – discrete frequency distribution

  • M (int) – number of resamplings

Return type:

1D array (integer)

Returns:

the resampled frequency distribution.

resample_diagline_dist(M)[source]

Return resampled frequency distribution of diagonal lines.

The resampled frequency distribution can be used for obtaining confidence bounds on diagonal line based RQA measures. This is described in detail in [Schinkel2009].

Concerning the choice of the number of resamplings, Schinkel et al. write: “The number of resamplings is not generally agreed upon but common guidelines suggest values between 800 and 1500.”

Parameters:

M (int) – number of resamplings

Return type:

1D array (integer)

Returns:

the resampled frequency distribution of diagonal lines.

resample_vertline_dist(M)[source]

Return resampled frequency distribution of vertical lines.

The resampled frequency distribution can be used for obtaining confidence bounds on vertical line based RQA measures. This is described in detail in [Schinkel2009].

Concerning the choice of the number of resamplings, Schinkel et al. write: “The number of resamplings is not generally agreed upon but common guidelines suggest values between 800 and 1500.”

Parameters:

M (int) – number of resamplings

Return type:

1D array (integer)

Returns:

the resampled frequency distribution of vertical lines.

rqa_summary(l_min=2, v_min=2)[source]

Return a selection of RQA measures.

The selection consists of the recurrence rate \(RR\), the determinism \(DET\), the average diagonal line length \(L\) and the laminarity \(LAM\).

Parameters:
  • l_min (int) – The minimum diagonal line length.

  • v_min (int) – The minimum vertical line length.

Return type:

Python dictionary

Returns:

a selection of RQA measures.

set_adaptive_neighborhood_size(adaptive_neighborhood_size, order=None)[source]

Construct recurrence plot using the adaptive neighborhood size algorithm introduced in [Xu2008].

The exact algorithm was deduced from private correspondence with the authors, as the description given in the above mentioned is not correct or at least ambiguous.

Modifies / sets the class variables R and N accordingly.

Parameters:
  • adaptive_neighborhood_size (number) – The number of adaptive nearest neighbors (recurrences) assigned to each state vector.

  • order (1D array of int32) – The indices of state vectors in the order desired for processing by the algorithm. The standard order is \(1,...,N\).

set_fixed_local_recurrence_rate(local_recurrence_rate)[source]

Set the recurrence plot to a fixed local recurrence rate.

This results in a fixed number of recurrences for each state vector, i.e., all state vectors have the same number of recurrences. Modifies / sets the class variables R and N accordingly.

Note

The resulting recurrence matrix \(R\) is generally asymmetric!

Parameters:

local_recurrence_rate (number) – The local recurrence rate.

set_fixed_recurrence_rate(recurrence_rate)[source]

Set the recurrence plot to a fixed recurrence rate.

Modifies / sets the class variables R and N accordingly.

Parameters:

recurrence_rate (number) – The recurrence rate.

set_fixed_threshold(threshold)[source]

Set the recurrence plot to a fixed threshold.

Modifies / sets the class variables R and N accordingly.

Parameters:

threshold (number) – The recurrence threshold.

set_fixed_threshold_std(threshold_std)[source]

Set the recurrence plot to a fixed threshold in units of the standard deviation of the time series.

Calculates the absolute threshold and calls set_fixed_threshold().

Parameters:

threshold_std (number) – The recurrence threshold in units of the standard deviation of the time series.

silence_level

The inverse level of verbosity of the object.

sparse_rqa

Controls sequential calculation of RQA measures.

supremum_distance_matrix()[source]

Return the supremum distance matrix from an embedding of a time series.

Return type:

2D square Numpy array

Returns:

the supremum distance matrix.

static threshold_from_recurrence_rate(distance, recurrence_rate: float)[source]

Return the threshold for recurrence plot construction given the recurrence rate.

Be aware, that the returned threshold can only approximately give the desired recurrence rate. The accuracy depends on the distribution of values in the given distance matrix \(D\).

Parameters:
  • distance (2D square array.) – The phase space distance matrix \(D\).

  • recurrence_rate (number) – The desired recurrence rate.

Return number:

the recurrence threshold corresponding to the desired recurrence rate.

static threshold_from_recurrence_rate_fast(distance, recurrence_rate, rr_precision=0.001)[source]

Return the threshold for recurrence plot construction given the recurrence rate.

The threshold yielding a given recurrence_rate is approximated using a randomly selected rr_precision percent of the distance matrix’ entries. Hence, the expected accuracy is lower than that achieved by using threshold_from_recurrence_rate.

Parameters:
  • distance (2D square array.) – The phase space distance matrix \(D\).

  • recurrence_rate (number) – The desired recurrence rate.

  • rr_precision (number) – The desired precision of recurrence rate estimation.

Return number:

the recurrence threshold corresponding to the desired recurrence rate.

time_series

The time series from which the recurrence plot is constructed.

trapping_time(v_min=2, resampled_dist=None)[source]

Alias for average_vertlength() (see description there).

twin_surrogates(n_surrogates=1, min_dist=7)[source]

Generate surrogates based on the current (embedded) time series embedding using the twin surrogate method.

The twins surrogates have the same dimensionality as the (embedded) trajectory used for constructing the recurrence plot. If scalar surrogate time series are desired, any component of the twin surrogate trajectory may be isolated.

Twin surrogates share linear and nonlinear properties with the original time series, since they correspond to realizations of trajectories of the same dynamical systems with different initial conditions.

References: [Thiel2006] [*], [Marwan2007].

Parameters:
  • min_dist (number) – The minimum temporal distance for twins.

  • n_surrogates (int) – The number of twin surrogate trajectories to be returned.

Return type:

3D array (surrogate number, time, dimension)

Returns:

the twin surrogate trajectories.

twins(min_dist=7)[source]

Return list of the twins of each state vector based on the recurrence matrix.

Two state vectors are said to be twins if they share the same recurrences, i.e., if the corresponding rows or columns in the recurrence plot are identical.

References: [Thiel2006], [Marwan2007].

Parameters:

min_dist (number) – The minimum temporal distance for twins.

Return [[number]]:

the list of twins for each state vector in the time series.

vert_entropy(v_min=2, resampled_dist=None)[source]

Return vertical line-based RQA measure vertical line entropy.

It is defined as the entropy of the probability to find a vertical line of exactly length l in the RP - reflects the complexity of the RP with respect to vertical lines.

Parameters:
  • v_min (int) – The minimal vertical line length.

  • resampled_dist (1D array (integer)) – resampled frequency distribution of vertical lines

Return number:

the vertical line-based entropy.

vertline_dist()[source]

Return the frequency distribution of vertical line lengths \(P(v-1)\).

Note that entry \(P(v-1)\) contains the number of vertical lines of length \(v\). Thus, \(P(0)\) counts lines of length \(1\), \(P(1)\) counts lines of length \(2\), asf.

Return type:

1D array (int32)

Returns:

the frequency distribution of vertical line lengths \(P(v-1)\).

white_vert_entropy(w_min=1)[source]

Return white vertical line-based RQA measure white vertical line entropy.

It is defined as the entropy of the probability to find a white vertical line of exactly length l in the RP - reflects the complexity of the RP with respect to white vertical lines (recurrence times).

Parameters:

w_min (int) – Minimal white vertical line length (recurrence time).

Return number:

the white vertical line-based entropy.

white_vertline_dist()[source]

Return the frequency distribution of white vertical line lengths \(P(w-1)\).

Note that entry \(P(w-1)\) contains the number of white vertical lines of length \(w\). Thus, \(P(0)\) counts lines of length \(1\), \(P(1)\) counts lines of length \(2\), asf.

The length of a white vertical line in a recurrence plot corresponds to the time the system takes to return close to an earlier state.

Return type:

1D array (int32)

Returns:

the frequency distribution of white vertical line lengths \(P(w-1)\).