API

Input data for RecurrenceMicrostatesAnalysis.jl

RecurrenceMicrostatesAnalysis.jl accepts three types of input, each associated with a different backend:

  • StateSpaceSet or Vector{<:Real}: used for multivariate time series, datasets, or state-space representations.
  • AbstractArray{<:Real}: used for spatial data, enabling RMA to be applied within the generalized framework of Spatial Recurrence Plots (SRP) (Marwan et al., 2007). We give some examples about its use in Spatial data.
  • AbstractGPUVector: used for time series analysis with the GPU backend. We provide some explanations about it in GPU.
Spatial Recurrence Microstates Analysis

RMA with SRP is an open research field. We include this functionality in the package for exploratory purposes, but the method is not yet mature enough for production use. Nevertheless, feel free to experiment with it in your research. 😃

StateSpaceSets.StateSpaceSetType
StateSpaceSet{D, T, V} <: AbstractVector{V}

A dedicated interface for sets in a state space. It is an ordered container of equally-sized points of length D, with element type T, represented by a vector of type V. Typically V is SVector{D,T} or Vector{T} and the data are always stored internally as Vector{V}. SSSet is an alias for StateSpaceSet.

The underlying Vector{V} can be obtained by vec(ssset), although this is almost never necessary because StateSpaceSet subtypes AbstractVector and extends its interface. StateSpaceSet also supports almost all sensible vector operations like append!, push!, hcat, eachrow, among others. When iterated over, it iterates over its contained points.

Construction

Constructing a StateSpaceSet is done in three ways:

  1. By giving in each individual columns of the state space set as Vector{<:Real}: StateSpaceSet(x, y, z, ...).
  2. By giving in a matrix whose rows are the state space points: StateSpaceSet(m).
  3. By giving in directly a vector of vectors (state space points): StateSpaceSet(v_of_v).

All constructors allow for two keywords:

  • container which sets the type of V (the type of inner vectors). At the moment options are only SVector, MVector, or Vector, and by default SVector is used.
  • names which can be an iterable of length D whose elements are Symbols. This allows assigning a name to each dimension and accessing the dimension by name, see below. names is nothing if not given. Use StateSpaceSet(s; names) to add names to an existing set s.

Description of indexing

When indexed with 1 index, StateSpaceSet behaves exactly like its encapsulated vector. i.e., a vector of vectors (state space points). When indexed with 2 indices it behaves like a matrix where each row is a point.

In the following let i, j be integers, typeof(X) <: AbstractStateSpaceSet and v1, v2 be <: AbstractVector{Int} (v1, v2 could also be ranges, and for performance benefits make v2 an SVector{Int}).

  • X[i] == X[i, :] gives the ith point (returns an SVector)
  • X[v1] == X[v1, :], returns a StateSpaceSet with the points in those indices.
  • X[:, j] gives the jth variable timeseries (or collection), as Vector
  • X[v1, v2], X[:, v2] returns a StateSpaceSet with the appropriate entries (first indices being "time"/point index, while second being variables)
  • X[i, j] value of the jth variable, at the ith timepoint

In all examples above, j can also be a Symbol, provided that names has been given when creating the state space set. This allows accessing a dimension by name. This is provided as a convenience and it is not an optimized operation, hence recommended to be used primarily with X[:, j::Symbol].

Use Matrix(ssset) or StateSpaceSet(matrix) to convert. It is assumed that each column of the matrix is one variable. If you have various timeseries vectors x, y, z, ... pass them like StateSpaceSet(x, y, z, ...). You can use columns(dataset) to obtain the reverse, i.e. all columns of the dataset in a tuple.

source

Output data from RecurrenceMicrostatesAnalysis.jl

When computing the RMA distribution, RecurrenceMicrostatesAnalysis.jl returns a Probabilities or Counts. This type is provided by ComplexityMeasures.jl, allowing this package to interoperate naturally with its tools and workflows.

ComplexityMeasures.ProbabilitiesType
Probabilities <: Array{<:AbstractFloat, N}
Probabilities(probs::Array [, outcomes [, dimlabels]]) → p
Probabilities(counts::Counts [, outcomes [, dimlabels]]) → p

Probabilities stores an N-dimensional array of probabilities, while ensuring that the array sums to 1 (normalized probability mass). In most cases the array is a standard vector. p itself can be manipulated and iterated over, just like its stored array.

The probabilities correspond to outcomes that describe the axes of the array. If p isa Probabilities, then p.outcomes[i] is an an abstract vector containing the outcomes along the i-th dimension. The outcomes have the same ordering as the probabilities, so that p[i][j] is the probability for outcome p.outcomes[i][j]. The dimensions of the array are named, and can be accessed by p.dimlabels, where p.dimlabels[i] is the label of the i-th dimension. Both outcomes and dimlabels are assigned automatically if not given. If the input is a set of Counts, and outcomes and dimlabels are not given, then the labels and outcomes are inherited from the counts.

Examples

julia> probs = [0.2, 0.2, 0.2, 0.2]; Probabilities(probs) # will be normalized to sum to 1
 Probabilities{Float64,1} over 4 outcomes
 Outcome(1)  0.25
 Outcome(2)  0.25
 Outcome(3)  0.25
 Outcome(4)  0.25
julia> c = Counts([12, 16, 12], ["out1", "out2", "out3"]); Probabilities(c)
 Probabilities{Float64,1} over 3 outcomes
 "out1"  0.3
 "out2"  0.4
 "out3"  0.3
source
ComplexityMeasures.CountsType
Counts <: Array{<:Integer, N}
Counts(counts [, outcomes [, dimlabels]]) → c

Counts stores an N-dimensional array of integer counts corresponding to a set of outcomes. This is typically called a "frequency table" or "contingency table".

If c isa Counts, then c.outcomes[i] is an abstract vector containing the outcomes along the i-th dimension, where c[i][j] is the count corresponding to the outcome c.outcomes[i][j], and c.dimlabels[i] is the label of the i-th dimension. Both labels and outcomes are assigned automatically if not given. c itself can be manipulated and iterated over like its stored array.

source

Recurrence Microstates

RecurrenceMicrostatesAnalysis.RecurrenceMicrostatesType
RecurrenceMicrostates <: CountBasedOutcomeSpace

An OutcomeSpace from ComplexityMeasures.jl representing recurrence microstates.

Description

Let $\mathscr{X} = (\vec{x}_i)_{i = 1}^{K}$ be a sequence of data elements. From $\mathscr{X}$, we can extract subsequences $\mathbf{X}_{(p)}$ of length $N$ such that $\mathbf{X}_{(p)} = (\vec{x}_i)_{i=p+1}^{p+N}$ for $0 \leq p \leq K - N$.

Let $r_{(i,j)}$ be a recurrence function that takes two elements, $\vec{x}_i$ and $\vec{x}_j$, and returns whether they are recurrent. For example, the threshold-based recurrence is given by $r_{(i,j)} = \Theta(\varepsilon - \|\vec{x}_i - \vec{x}_j\|)$, where $\Theta$ is the Heaviside function, $\|\cdot\|$ is an appropriate metric for the data, and $\varepsilon$ is the recurrence threshold. In this representation, a recurrence is denoted by 1 and a non-recurrence by 0.

For two subsequences belonging to $\mathscr{X}$, we can generate, using the recurrence function, a matrix $\mathbf{M} \equiv \mathbf{R}(\mathbf{X}_{(p)}, \mathbf{X}_{(q)})$:

\[\mathbf{R}(\mathbf{X}_{(p)}, \mathbf{X}_{(q)}) \coloneqq \begin{bmatrix} r_{p+1,q+1} & \dots & r_{p+1,q+N} \\ \vdots & \ddots & \vdots \\ r_{p+N,q+1} & \dots & r_{p+N,q+N} \end{bmatrix}.\]

This matrix is known as a recurrence microstate (Corso et al., 2018) of length $N$ when $N \ll K$, and represents a local portion of the recurrence plot $\mathbf{R}(\mathscr{X}, \mathscr{X})$.

Implementation

To define a recurrence microstate, three components are required:

  1. The recurrence function, $r_{(i,j)}$, used to compute recurrences (see RecurrenceExpression).
  2. The microstate shape and its size (see MicrostateShape).
  3. The method used to extract microstates from $\mathscr{X}$ (see SamplingMode).

Constructors

RecurrenceMicrostates(expr::RecurrenceExpression, shape::MicrostateShape; kwargs...)
RecurrenceMicrostates(expr::RecurrenceExpression, N::Int; kwargs...) # It uses square microstates.

Using ThresholdRecurrence:

RecurrenceMicrostates(ε::Real, N::Int; kwargs...)
RecurrenceMicrostates(ε::Real, shape::MicrostateShape; kwargs...)

Using CorridorRecurrence:

RecurrenceMicrostates(ε_min::Real, ε_max::Real, N::Int; kwargs...)
RecurrenceMicrostates(ε_min::Real, ε_max::Real, shape::MicrostateShape; kwargs...)

Spatial generalization

We also implement a spatial generalization of recurrence microstate analysis based on (Marwan et al., 2007). It operates similarly to standard RMA, but retrieves microstates from a recurrence tensor constructed from the data (without explicitly constructing the full tensor). This is an experimental feature included as an exploratory addition to the package; there is currently no established literature describing how to use this generalization. If you are curious, feel free to experiment with it 🙂

In this setting, the input data is no longer a sequence of elements, but a spatial structure with $d$ dimensions. Each element of this structure can be either a Real or a Vector containing features associated with the corresponding position. Following (Marwan et al., 2007), the resulting recurrence space has $2 \times d$ dimensions. For example, for an image we have $d = 2$, so the recurrence space has 4 dimensions.

The recurrence expression is reformulated to operate on vector indices: $r_{(\vec{i}, \vec{j})}$. For example, the threshold-based recurrence is given by $r_{(\vec{i}, \vec{j})} = \Theta(\varepsilon - \|\mathbf{x}_{\vec{i}} - \mathbf{x}_{\vec{j}}\|)$, where the coordinate vectors are defined as:

\[\vec{i} = \sum_{q = 1}^{d} i_q \hat{e}_q,\]

where $\hat{e}_q$ is the unit vector representing the corresponding dimension $q$.

Using this formulation, recurrence microstates are defined as before. However, since the recurrence plot is now a tensor, microstates may have hypergeometric shapes, such as hypercubes or hyperrectangles, or correspond to projections onto lower-dimensional subspaces. This provides greater flexibility in capturing structure, but also increases the complexity of using the method.

Moreover, since RQA is defined differently for spatial recurrence plots, the implementations in RecurrenceMicrostatesAnalysis.jl are generally not compatible with this generalization, except for RecurrenceRate and Entropy, which can be estimated from any recurrence microstate distribution.

Shape and sampling compatibility

Some shapes and sampling modes are not compatible with the spatial generalization, e.g. TriangleMicrostate and Full.

Spatial constructors

These constructors use a hypergeometric version of RectMicrostate, defined by an NTuple called structure.

RecurrenceMicrostates(expr::RecurrenceExpression, structure::NTuple; kwargs...)
RecurrenceMicrostates(ε::Real, structure::NTuple; kwargs...)
RecurrenceMicrostates(ε_min::Real, ε_max::Real, structure::NTuple; kwargs...)

Keyword arguments

  • metric: The metric used to compute distances between points. Typically it is a subtype of Metric but it can be an arbitrary function of two points returning a real. In the case of using GPUs, this must be an instance of GPUMetric.
  • ratio: The sampling ratio. The default is 0.1.
  • sampling: The sampling mode. The default is SRandom.
Memory

Note that the number of possible microstates is $2^\sigma$, where $\sigma$ is the number of recurrence entries in the microstate structure. This means that as the number of entries increases, the memory required to store the full distribution quickly becomes impractical. For example, a square $6 \times 6$ microstate has 36 entries, resulting in $2^{36} = 68,719,476,736$ possible microstates.

source

Recurrence expressions

RecurrenceMicrostatesAnalysis.recurrenceFunction
recurrence(expr::RecurrenceExpression, [x], [y], [...]) → UInt8

Defines how the recurrence state between [x] and [y] is computed for a given RecurrenceExpression.

The additional arguments ([...]) depend on whether recurrence is computed for time series or spatial data.

Time series recurrence

function recurrence(expr::RecurrenceExpression, x::StateSpaceSet, y::StateSpaceSet, ::Int, ::Int)

The two Int arguments correspond to the positions $(i, j)$ in the time series used to evaluate recurrence.

Spatial recurrence

function recurrence(expr::RecurrenceExpression, x::AbstractArray{<:Real}, y::AbstractArray{<:Real}, ::NTuple{N, Int}, ::NTuple{M, Int})

The two NTuple{N, Int} arguments correspond to the positions $(\vec{i}, \vec{j})$ in the spatial data used to evaluate recurrence.

Info

To support GPU execution, recurrence expressions must implement gpu_recurrence instead of recurrence. The arguments are equivalent, with the addition of the phase-space dimension as an input parameter. See ThresholdRecurrence for a reference implementation.

Output

The recurrence function must always return a UInt8: 0 for non-recurrence and 1 for recurrence.

source
RecurrenceMicrostatesAnalysis.ThresholdRecurrenceType
ThresholdRecurrence <: RecurrenceExpression

Recurrence expression defined by the standard threshold criterion:

\[r_{(i,j)} = \Theta(\varepsilon - |\vec{x}_i - \vec{y}_j|),\]

where $\Theta(\cdot)$ denotes the Heaviside function and $\varepsilon$ is the distance threshold defining the maximum separation for two states to be considered recurrent.

The ThresholdRecurrence struct stores the threshold parameter ε, as well as the distance metric used to evaluate $\|\vec{x} - \vec{y}\|$. The metric must be defined using the Distances.jl package.

If the data for $\vec{x}$ and $\vec{y}$ are the same, the result is a recurrence plot; otherwise, it is a cross-recurrence plot.

Constructor

ThresholdRecurrence(ε::Real; metric::Metric = Euclidean())

Examples

ThresholdRecurrence(0.27)
ThresholdRecurrence(0.27; metric = Cityblock())
ThresholdRecurrence(0.27; metric = GPUEuclidean())

The recurrence evaluation is performed via the recurrence function. For GPU execution, the corresponding implementation is provided by gpu_recurrence.

source
RecurrenceMicrostatesAnalysis.CorridorRecurrenceType
CorridorRecurrence <: RecurrenceExpression

Recurrence expression defined by the corridor criterion introduced in (Iwanski and Bradley, 1998):

\[r_{(i, j)} = \Theta(|\vec{x}_i - \vec{y}_j| - \varepsilon_{min}) \cdot \Theta(\varepsilon_{max} - |\vec{x}_i - \vec{y}_j|),\]

where $\Theta(\cdot)$ denotes the Heaviside function and $(\varepsilon_{\min}, \varepsilon_{\max})$ define the minimum and maximum distance thresholds for two states to be considered recurrent.

The CorridorRecurrence struct stores the corridor thresholds ε_min and ε_max, as well as the distance metric used to evaluate $\|\vec{x}_i - \vec{y}_j\|$. The metric must be defined using the Distances.jl package.

If the data for $\vec{x}$ and $\vec{y}$ are the same, the result is a recurrence plot; otherwise, it is a cross-recurrence plot.

Constructor

CorridorRecurrence(ε_min::Real, ε_max::Real; metric::Metric = Euclidean())

Examples

CorridorRecurrence(0.05, 0.27)
CorridorRecurrence(0.05, 0.27; metric = Cityblock())
CorridorRecurrence(0.05, 0.27; metric = GPUEuclidean())

The recurrence evaluation is performed via the recurrence function. For GPU execution, the corresponding implementation is provided by gpu_recurrence.

source

Microstate shapes

RecurrenceMicrostatesAnalysis.MicrostateShapeType
MicrostateShape

Abstract supertype defining the basic structure and layout of a microstate (or motif).

A MicrostateShape specifies which relative positions are retrieved from the data to evaluate recurrences, and how these binary recurrence values are interpreted and mapped to a decimal representation for counting.

All subtypes of MicrostateShape must include a field expr, which defines the RecurrenceExpression used to compute recurrences.

Implementations

source
RecurrenceMicrostatesAnalysis.RectMicrostateType
RectMicrostate <: MicrostateShape

Define a rectangular microstate shape.

RectMicrostate can represent either a two-dimensional microstate (identified as Rect2Microstate, used for Recurrence Plots and Cross-Recurrence Plots) or an N-dimensional microstate (identified as RectNMicrostate, used for spatial data).

Rect2Microstate (time-series data)

A 2D rectangular microstate can be initialized using either of the following constructor:

RectMicrostate(rows::Int, cols::Int; B = 2)

Here, rows and columns define the rectangle dimensions, and B is the base used to encode microstate elements (typically 2, representing recurrence or non-recurrence).

Rectangular microstates can be specialized to define common patterns such as lines, columns, and squares:

line = RectMicrostate(N, 1)
column = RectMicrostate(1, N)
square = RectMicrostate(N, N)

Since square microstates are frequently used, a convenience constructor is also provided:

RectMicrostate(N; B = 2)

RectNMicrostate (spatial data)

For N-dimensional structures, typically used with spatial data, the RectNMicrostate variant can be initialized as:

RectMicrostate(structure::NTuple{D, Int}; B = 2)

Here, structure defines the size of the microstate along each dimension. For example:

nrect = RectMicrostate((2, 1, 2, 1))

This form is suitable for N-dimensional spatial data, such as images or volumetric datasets.

source
RecurrenceMicrostatesAnalysis.DiagonalMicrostateType
DiagonalMicrostate <: MicrostateShape

Define a diagonal microstate shape, which captures recurrences along diagonals of a Recurrence Plot (RP).

Constructor

DiagonalMicrostate(N::Int; B::Int = 2)

where N defines the length of the diagonal microstate.

Example

diagonal = DiagonalMicrostate(3)
Info

DiagonalMicrostate microstates are compatible with spatial data. However, they do not capture hyper-diagonals in Spatial Recurrence Plots (SRP). Only diagonals defined by sequential recurrences are supported, such as:

\[R_{i_1,i_2,j_1,j_2}, R_{i_1 + 1,i_2 + 1,j_1 + 1,j_2 + 1}, R_{i_1 + 2,i_2 + 2,j_1 + 2,j_2 + 2}, \ldots, R_{i_1 + n - 1,i_2 + n - 1,j_1 + n - 1,j_2 + n - 1}\]

source
RecurrenceMicrostatesAnalysis.TriangleMicrostateType
TriangleMicrostate <: MicrostateShape

Define a triangular microstate shape, originally introduced by Hirata in 2021 (Hirata, 2021).

Constructor

TriangleMicrostate(N::Int; B::Int = 2)

where N defines the size of the triangular microstate.

Example

N = 3
triangle = TriangleMicrostate(N)
Compat

Triangular microstate shape is not compatible with spatial data.

source

Sampling modes

RecurrenceMicrostatesAnalysis.SRandomType
SRandom{F<:Real} <: SamplingMode

Sampling mode that randomly selects microstate positions $(i, j)$ within the SamplingSpace.

Constructors

SRandom(num_samples::Int)
SRandom(ratio::Union{Float32, Float64})

The sampling mode can be initialized either by specifying the exact number of microstates to sample or by providing a ratio of the total number of possible microstates.

Examples

s = SRandom(1000)   # Specify the exact number of sampled microstates
s = SRandom(0.05)   # Specify a ratio of the total possible microstates
source

Sampling space

Future implementation

We pretend to expand the RecurrenceMicrostates structure to also consider a setted space from the recurrence plot as source of information to construct the RMA distribution.

Recurrence Quantification Analysis

RecurrenceMicrostatesAnalysis.RecurrenceRateType
RecurrenceRate <: ComplexityEstimator
RecurrenceRate(ε::Float64, N::Int = 1; kwargs...)

An estimator of the recurrence rate, used with complexity. It uses a $1 \times 1$ microstate by default, but you can set a different size via the N parameter. The recurrence rate is estimated for a threshold ε.

Keyword arguments

  • metric::Metric: The metric used to compute recurrence.
  • ratio: The sampling ratio. The default is 0.1.
  • sampling: The sampling mode. The default is SRandom.

Description

Recurrence rate (RR) is defined as (Webber and Marwan, 2015)

\[RR = \frac{1}{K^2}\sum_{i,j=1}^{K} r_{(i,j)}.\]

When estimating it using RMA, the recurrence rate is defined as the expected value over the microstate distribution:

\[RR \approx \sum_{i=1}^{2^\sigma} p_i^{(N)}RR_i^{(N)},\]

where $RR_i^{(N)}$ denotes the recurrence rate of the $i$-th microstate.

Performance

Although estimating RR using RMA is faster than typical RQA computation, the precision depends on the time series length. Therefore, for small time series, i.e., $K \leq 1000$, we strongly recommend using standard RQA with RecurrenceAnalysis.jl.

source
RecurrenceMicrostatesAnalysis.RecurrenceDeterminismType
RecurrenceDeterminism <: ComplexityEstimator
RecurrenceDeterminism(ε::Float64; kwargs...)

An estimator of recurrence determinism, used with complexity. Determinism is estimated for a threshold ε.

Keyword arguments

  • metric::Metric: The metric used to compute recurrence.
  • ratio: The sampling ratio. The default is 0.1.
  • sampling: The sampling mode. The default is SRandom.

Description

Recurrence determinism (DET) is defined as (Webber and Marwan, 2015)

\[DET = \frac{\sum_{l=l_{min}}^{K} l H_D(l)}{\sum_{i,j=1}^{K} r_{(i,j)}},\]

where $H_D(l)$ is the histogram of diagonal line lengths:

\[H_D(l) = \sum_{i,j=1}^{K} (1 - r_{i-1, j-1})(1 - r_{i+l,j+l})\prod_{k=0}^{l-1} r_{i+k,j+k}.\]

By inverting the determinism expression, we can rewrite it as (da Cruz et al., 2025)

\[DET = 1 - \frac{1}{K^2 \sum_{i,j=1}^{K} r_{(i,j)}} \sum_{l=1}^{l_{min} - 1} l H_D(l).\]

An approximate value for DET can be estimated using recurrence microstates, as introduced by da Cruz et al. (da Cruz et al., 2025). From an input dataset x, we estimate a recurrence microstate distribution $\vec{p}$. This distribution must be defined over square microstates of size $3 \times 3$. Here, we use the relation:

\[\frac{H_D(l)}{(K-l-1)^2} = \vec{d}^{(l)} \cdot \mathcal{R}^{(l + 2)}\vec{p}^{(l + 2)}.\]

For the commonly used case $l_{min} = 2$, this leads to the approximation

\[DET \approx 1 - \frac{\vec{d}^{(1)}\cdot\mathcal{R}^{(3)}\vec{p}^{(3)}}{\sum_{i,j=1}^{K} r_{(i,j)}}.\]

The correlation term $\vec{d}^{(1)} \cdot \mathcal{R}^{(3)} \vec{p}^{(3)}$ can be simplified by explicitly identifying the microstates selected by $\vec{d}^{(1)}$. This corresponds to selecting only microstates of the form:

\[\begin{pmatrix} \xi & \xi & 0 \\ \xi & 1 & \xi \\ 0 & \xi & \xi \end{pmatrix},\]

where $\xi$ denotes an unconstrained entry. There are 64 microstates with this structure among the 512 possible $3 \times 3$ microstates. Defining the class $C_D$ as the set of microstates with this structure, DET can be estimated as:

\[DET \approx 1 - \frac{\sum_{i\in C_D} p_i^{(3)}}{\sum_{i,j=1}^{K} r_{(i,j)}}.\]

The implementation used by complexity is an optimized version of this process using DiagonalMicrostate (Ferreira et al., 2025). Since this microstate shape is symmetric with respect to the desired information, it is not necessary to account for $\xi$ values as in the square microstate case. Thus, determinism can be estimated as

\[DET \approx 1 - \frac{p_3^{(3)}}{\sum_{i,j=1}^{K} r_{(i,j)}},\]

where $p_3^{(3)}$ is the probability of observing the microstate $0~1~0$.

Performance

Although estimating DET using RMA is faster than typical RQA computation, the precision depends on the time series length. Therefore, for small time series, i.e., $K \leq 1000$, we strongly recommend using standard RQA with RecurrenceAnalysis.jl.

source
RecurrenceMicrostatesAnalysis.RecurrenceLaminarityType
RecurrenceLaminarity <: ComplexityEstimator
RecurrenceLaminarity(ε::Float64; kwargs...)

An estimator of recurrence laminarity, used with complexity. Laminarity is estimated for a threshold ε.

Keyword arguments

  • metric::Metric: The metric used to compute recurrence.
  • ratio: The sampling ratio. The default is 0.1.
  • sampling: The sampling mode. The default is SRandom.

Description

Recurrence laminarity (LAM) is defined as (Webber and Marwan, 2015)

\[LAM = \frac{\sum_{l=l_{min}}^{K} l H_V(l)}{\sum_{i,j=1}^{K} r_{(i,j)}},\]

where $H_V(l)$ is the histogram of vertical line lengths:

\[H_V(l) = \sum_{i,j=1}^{K} (1 - r_{i, j-1})(1 - r_{i,j+l})\prod_{k=0}^{l-1} r_{i,j+k}.\]

By inverting the laminarity expression, we can rewrite it as (da Cruz et al., 2025)

\[LAM = 1 - \frac{1}{K^2 \sum_{i,j=1}^{K} r_{(i,j)}} \sum_{l=1}^{l_{min} - 1} l H_V(l).\]

An approximate value for LAM can be estimated using recurrence microstates, as introduced by da Cruz et al. (da Cruz et al., 2025). From an input dataset x, we estimate a recurrence microstate distribution $\vec{p}$. This distribution must be defined over square microstates of size $3 \times 3$. Here, we use the relation:

\[\frac{H_V(l)}{(K-l-1)^2} = \vec{v}^{(l)} \cdot \mathcal{R}^{(l + 2)}\vec{p}^{(l + 2)}.\]

For the commonly used case $l_{min} = 2$, this leads to the approximation

\[LAM \approx 1 - \frac{\vec{v}^{(1)}\cdot\mathcal{R}^{(3)}\vec{p}^{(3)}}{\sum_{i,j=1}^{K} r_{(i,j)}}.\]

The correlation term $\vec{v}^{(1)} \cdot \mathcal{R}^{(3)} \vec{p}^{(3)}$ can be simplified by explicitly identifying the microstates selected by $\vec{v}^{(1)}$. This corresponds to selecting only microstates of the form:

\[\begin{pmatrix} 0 & 1 & 0 \\ \xi & \xi & \xi \\ \xi & \xi & \xi \end{pmatrix},\]

where $\xi$ denotes an unconstrained entry. There are 64 microstates with this structure among the 512 possible $3 \times 3$ microstates. Defining the class $C_V$ as the set of microstates with this structure, LAM can be estimated as:

\[LAM \approx 1 - \frac{\sum_{i\in C_V} p_i^{(3)}}{\sum_{i,j=1}^{K} r_{(i,j)}}.\]

The implementation used by complexity is an optimized version of this process using $1 \times 3$ RectMicrostate (Ferreira et al., 2025). Since this microstate shape is symmetric with respect to the desired information, it is not necessary to account for $\xi$ values as in the square microstate case. Thus, laminarity can be estimated as

\[LAM \approx 1 - \frac{p_3^{(3)}}{\sum_{i,j=1}^{K} r_{(i,j)}},\]

where $p_3^{(3)}$ is the probability of observing the microstate $0~1~0$.

Performance

Although estimating LAM using RMA is faster than typical RQA computation, the precision depends on the time series length. Therefore, for small time series, i.e., $K \leq 1000$, we strongly recommend using standard RQA with RecurrenceAnalysis.jl.

source
RecurrenceMicrostatesAnalysis.DisorderType
Disorder{N} <: ComplexityEstimator
Disorder(N::Int = 4; kwargs...)

An estimator of a disorder measure, introduced by (Flauzino et al., 2025), used with complexity. It uses $N \times N$ microstates and a specified metric to compute the disorder.

Keyword arguments

  • metric::Metric: The metric used to compute recurrence.
  • threshold_range::Int: The number of threshold values used to estimate the disorder.

Description

Disorder, or the disorder index via symmetry in recurrence microstates (DISREM), is based on the implications of the disorder condition in microstates. It quantifies disorder through the information entropy of classes of recurrence microstates, which are required to be equiprobable within the same class according to the disorder condition.

Let $\sigma$ be a permutation, and $S_N$ the set of all permutations of $N$ elements, with $\sigma \in S_N$. Also, let $\mathcal{L}_\sigma$ be the operator that permutes the rows of a matrix, and $\mathcal{T}$ the operator that transposes a matrix. A recurrence microstate class is defined as (Flauzino et al., 2025):

\[\mathscr{M}_a (\mathbf{M}) = \bigcup_{\sigma_i,\sigma_j \in S_N} \{ \mathcal{L}_{\sigma_j}\mathcal{T}\mathcal{L}_{\sigma_i}\mathbf{M},\quad\mathcal{T}\mathcal{L}_{\sigma_j}\mathcal{T}\mathcal{L}_{\sigma_i}\mathbf{M} \}.\]

Let $p(\mathbf{M})$ be the probability of the microstate $\mathbf{M}$. We renormalize this probability with respect to its class:

\[p^{(a)}(\mathbf{M}) = \frac{p(\mathbf{M})}{\sum_{\mathbf{M}' \in \mathscr{M}_a} p(\mathbf{M}')}.\]

The normalized information entropy associated with the probability distribution of microstates in the class $\mathscr{M}_a$ is then defined as (Flauzino et al., 2025):

\[\xi_a(\varepsilon) = \frac{1}{m_a} \sum_{\mathbf{M} \in \mathscr{M}_a} p^{(a)}(\mathbf{M}) \ln p^{(a)}(\mathbf{M}),\]

where $m_a$ is the number of microstates in the class $\mathscr{M}_a$.

By summing the entropy over all classes and normalizing by its maximum amplitude $A$, we obtain the total entropy across all classes:

\[\xi(\varepsilon) = \frac{1}{A}\sum_{a = 1}^{A} \xi_a(\varepsilon).\]

The disorder measure is then defined as the maximum value of $\xi(\varepsilon)$ (Flauzino et al., 2025):

\[\Xi = \max_{\varepsilon} \xi(\varepsilon).\]

Microstate size and shape

Disorder is only defined for square microstates, with computations available for $N \in \{2, 3, 4, 5\}$ due to computational limitations. In particular, computations for $N = 5$ are computationally expensive.

source
RecurrenceMicrostatesAnalysis.WindowedDisorderType
WindowedDisorder{N, W} <: ComplexityEstimator
WindowedDisorder(W::Int, N::Int = 4; kwargs...)

This estimator is equivalent to Disorder, but computes it by splitting the data into windows of length W, returning a vector of measured disorder values for each window.

Keyword arguments

  • metric: The metric used to compute recurrence.
  • threshold_range::Int: The number of threshold values used to estimate disorder.
  • step::Int: The step between windows. The default is W.
source
RecurrenceMicrostatesAnalysis.rmaFunction
rma(ε::Float64, [x]; kwargs...) → Dict{Symbol, Float64}

Calculate all RMA estimations for a threshold ε and a time series x. All values are estimated using $3\times 3$ square microstates and a Full sampling mode.

Return

The returned value constains the following entries, which can be retrieved as a dictionary (e.g. results[:RR], etc.):

All the parameters returned by rma are Float64 numbers.

Keyword arguments

  • metric::Metric: The metric used to compute recurrence.
source

Optimization

RecurrenceMicrostatesAnalysis.optimizeFunction
optimize(param::Parameter, args...)

Optimize a free Parameter using a specific complexity measure.

Performance

The optimize function may compute multiple distributions and can be computationally expensive. Avoid calling it inside performance-critical loops.

source
RecurrenceMicrostatesAnalysis.ThresholdType
Threshold <: Parameter

The threshold is a free parameter used to classify two states as recurrent or non-recurrent. This parameter can be optimized using the optimize function in combination with specific complexity measures, e.g., Recurrence Entropy or Disorder.

Use:

optimize(p::Threshold, q::Entropy, N::Int, [x]; kwargs...)
optimize(p::Threshold, q::Disorder{N}, [x]; kwargs...)

Arguments

  • n: The size of the square microstate used in the optimization.
  • [x]: The input data.

Keyword arguments

  • ratio: The sampling ratio. The default is 0.1.
  • sampling: The sampling mode. The default is SRandom.
  • th_max_range: The fraction of the maximum distance defining the upper bound of the threshold search range. The default is 0.5.
  • th_start: The initial value of the threshold search range. The default is 1e-6.
  • fraction: The interaction fraction controlling the refinement process. The default is 5.
  • metric::Metric: The metric used to compute recurrence.
source

Operations

RecurrenceMicrostatesAnalysis.PermuteRowsType
PermuteRows{R, C} <: Operation

Operation that permutes the rows of a microstate $\mathbf{M}$.

To initialize a PermuteRows operation, a rectangular microstate shape must be provided via a RectMicrostate structure:

PermuteRows(::Rect2Microstate{R, C, B})

Examples

PermuteRows(RectMicrostate(3, 3))     #   Microstate 3 x 3
PermuteRows(RectMicrostate(3, 1))     #   Microstate 3 x 1 (it is a column)

This operation is applied via the operate function:

operate(::PermuteRows, I::Int, σ::Vector{Int})

Arguments

  • op: A PermuteRows operation.
  • I: Decimal identifier of the microstate (1-based).
  • σ: Permutation of rows to be applied.

Returns

The resulting microstate binary identifier (1-based).

source
RecurrenceMicrostatesAnalysis.PermuteColumnsType
PermuteColumns{R, C} <: Operation

Operation that permutes the columns of a microstate $\mathbf{M}$.

To initialize a PermuteColumns operation, a rectangular microstate shape must be provided via a RectMicrostate structure:

PermuteColumns(::Rect2Microstate{R, C, B}; S::Vector{Vector{Int}} = collect(permutations(1:C))

Here, the keyword argument S defines the set $S_n$ of column permutations. The PermuteColumns struct precomputes the column permutations for each row of the microstate. These precomputed permutations can be accessed via the field Q.

Examples

PermuteColumns(RectMicrostate(3, 3))     #   Microstate 3 x 3
PermuteColumns(RectMicrostate(1, 3))     #   Microstate 1 x 3 (it is a line)

This operation is applied via the operate function:

operate(op::PermuteColumns, I::Int, Qi::Int)

Arguments

  • op: A PermuteColumns operation.
  • I: Decimal identifier of the microstate (1-based).
  • Qi: Index of the permutation in the set S.

Returns

The resulting microstate decimal identifier (1-based).

source
RecurrenceMicrostatesAnalysis.TransposeType
Transpose{R, C} <: Operation

Operation that transposes a microstate $\mathbf{M}$.

To initialize a Transpose operation, a rectangular microstate shape must be provided via a RectMicrostate structure:

Transpose(::Rect2Microstate{R, C, B})

Examples

Transpose(RectMicrostate(3, 3))     # 3 x 3 microstate

This operation is applied via the operate function:

operate(::Transpose, I::Int)

Arguments

  • op: A Transpose operation.
  • I: Decima identifier of the microstate (1-based).

Returns

The resulting microstate decimal identifier (1-based).

source

Utils

GPU Metrics

Distances.jl is not compatible with GPU execution; therefore, distance evaluations must be reimplemented for GPU usage. This is done using a GPUMetric.

RecurrenceMicrostatesAnalysis.GPUMetricType
GPUMetric <: Metric

Abstract supertype for metrics compatible with the GPU backend.

Metrics subtyping GPUMetric must implement the internal evaluation function gpu_evaluate, which is used during GPU-based computations.

Implementations

source