Public API
Distribution
Main.RecurrenceMicrostatesAnalysis.distribution — Function
Based on Recurrence Plot
distribution([x], [parameters], n::Int; kwargs...)Compute the distribution of recurrence microstates probabilities from the dataset x. The input parameters consists of the constant values used to calculate the recurrence between two points. n is an integer that represents the length of motifs side.
Input:
[x]: input dataset.[parameter]: set of parameters used to compute the recurrence microstate distribution.n: microstate size.
Output: this function can return a vector, an array or a dictionary based on the number of possible microstates and the setting of run_mode or sampling_mode.
Based on Cross-Recurrence Plot
distribution([x], [y], parameters, n::Int; kwords...)Compute the distribution of recurrence microstates probabilities from the datasets x and y. The input parameters consists of the constant values used to calculate the recurrence between two points. n is an integer that represents the length of motifs side.
Input:
[x]: input dataset.[y]: input dataset.[parameter]: set of parameters used to compute the recurrence microstate distribution.n: microstate size.
Output: this function can return a vector, an array or a dictionary based on the number of possible microstates and the setting of run_mode or sampling_mode.
Using DifferencialEquations.jl
distribution([solution], parameters, n::Int, vicinity::Union{Int, Float64}; kwords...)Compute the distribution of recurrence microstates probabilities from the solution of a differencial equation solved by the library DifferencialEquations.jl. The input parameters consists of the constant values used to calculate the recurrence between two points. n is an integer that represents the length of motifs side. vicinity is the time separation used to discretize a continuous problem.
Input:
[solution]: solution returned by the libraryDifferentialEquations.jl.[parameter]: set of parameters used to compute the recurrence microstate distribution.n: microstate size.σ: sampling parameter; it defines the time resolution of discretized data.
Specific kwargs:
transient: defines an interval of time that will be ignored, and taked as a transient.K: defines the maximum size of the result time series.
Output: this function can return a vector, an array or a dictionary based on the number of possible microstates and the setting of run_mode or sampling_mode.
Main
distribution([x], [y], parameters, [structure]; kwords...)Compute the distribution of recurrence microstates probabilities from the datasets x and y. The input parameters consists of the constant values used to calculate the recurrence between two points. Meanwhile, the input structure is a vector where each element represents a side of the motif.
Input:
[x]: input dataset.[y]: input dataset.[parameter]: set of parameters used to compute the recurrence microstate distribution.[structure]: microstate structure.
kwargs:
shape: microstate shape. Can be:square,:triangle,:pair,:diagonalor:line. (default:square)run_mode: define the output format. It can be:vectfor aVector{Float64}, or:dictfor aDict{Int, Float64}. If you are you sampling_mode:columnwise`sampling_mode: define how the library will take motifs in a RP. Can be:full,:random,:triangleup,:columnwiseor:columnwise_full. (default:random)num_samples: number of samples used to compute the distribution. Can be anIntvalue or aFloat64, which will be interpretad as a proportion of the total population of microstates in a RP. (This is not required for:fulland:columnwise_fullsampling modes)threads: set if library will use asyncronous jobs or not.metric: metric defined using the libraryDistances.jl.func: recurrence function.
Output: this function can return a vector, an array or a dictionary based on the number of possible microstates and the setting of run_mode or sampling_mode.
RQA
Main.RecurrenceMicrostatesAnalysis.rrate — Function
rrate([probs])Compute the approximated recurrence rate of a RP from a probability distribution of recurrence microstates. Here, we use a relation between the mean recurrence rate of each motif and the desired value. It can be written as
\[RR \approx \sum_{I = 0}^N \mathbf{p}_I^{(k)}\left(\frac{1}{k^2}\sum_{i=1}^k\sum_{j=1}^k \mathbf{M}_{ij}^{(I)}\right),\]
where $\mathbf{M}_{ij}^{(I)}$ is the motif structure.
Input:
[probs]: the vector of probabilities $\mathbf{p}^{(k)}$ computed usingdistribution(...).
Output: returns the recurrence rate as a Float64.
rrate([x], [parameters], n::Int; shape::Symbol, sampling_mode::Symbol, r::Float64})Compute the approximated recurrence rate of a Recurrence Plot from a data [x] using a probability distribution of recurrence microstates computed from it.
Input:
[x]: input data.[parameter]: set of parameters used to compute the recurrence microstate distribution.n: microstate size.shape(kwarg): shape of the used motifs.:squareby default, it can be::square, :triangle, :pair, :diagonal, :line.sampling_mode(kwarg): sampling mode used.:randomby default, it can be::square, :triangle, :pair, :diagonal, :line.r(kwarg): ratio of the total number of microstates to be sampled for the histogram. (defaultr = 0.05)
Output: returns the recurrence rate as a Float64.
Main.RecurrenceMicrostatesAnalysis.rentropy — Function
Main.RecurrenceMicrostatesAnalysis.determinism — Function
determinism(rr::Float64, [probs])Estimate the determinism from a distribution. If the distribution has 512 elements, the function will consider square motifs, computing determinism using
\[I^{(\beta)} = \frac{1}{RR}(2\mathbf{M}_{1,2}^{(\beta)} + 4\mathbf{M}_{1,3}^{(\beta)} + 8\mathbf{M}_{2,1}^{(\beta)} + 16 + 32\mathbf{M}_{2,3}^{(\beta)} + 64\mathbf{M}_{3,1}^{(\beta)} + 128\mathbf{M}_{3,2}^{(\beta)}),\]
where $\mathbf{M}$ is the motif structure. It defines a class of motifs $(C_L) \ni I^{(\beta)}$ that we use to estimate DET:
\[DET \approx 1 - \frac{1}{RR}\sum_{I\in (C_L)} \mathbf{p}^{(3)}_I.\]
dist = distribution(data, th, 3)
rr = rrate(dist)
det = determinism(rr, dist)If the distribution has 8 elements, this function will consider line motifs, which makes the process simpler. In this case, we just need the motif with $I = 2$:
\[DET \approx 1 - \frac{\mathbf{p}^{(3)}_2}{RR}.\]
dist = distribution(data, th, 3; shape = :line)
rr = rrate(dist)
det = determinism(rr, dist)Input:
rr: recurrence rate.[probs]: aVector{Float64}returned by the functiondistribution(...).
Output: returns the determinism as a Float64.
determinism([x], threshold::Float64; r::Float64)Estimate the determinism from a data [x]` using a probability distribution and a RR computed from it.
Input:
[x]: input data.[parameter]: set of parameters used to compute the recurrence microstate distribution.r(kwarg): ratio of the total number of microstates to be sampled for the histogram. (defaultr = 0.05)
Output: returns a Tuple{Float64, Float64}.
lam: laminarity asFloat64.rr: recurrence rate asFloat64.
Main.RecurrenceMicrostatesAnalysis.laminarity — Function
laminarity(rr::Float64, [probs])Estimate the laminarity from a distribution. If the distribution has 512 elements, the function will consider square motifs, computing laminarity using
\[I^{(\beta)} = \frac{1}{RR}(2 + 8\mathbf{M}_{2,1}^{(\beta)} + 16\mathbf{M}_{2,2}^{(\beta)} + 32\mathbf{M}_{2,3}^{(\beta)} + 64\mathbf{M}_{3,1}^{(\beta)} + 128\mathbf{M}_{3,2}^{(\beta)} + 256\mathbf{M}_{3,3}^{(\beta)}),\]
where $\mathbf{M}$ is the motif structure. It defines a class of motifs $(C_L) \ni I^{(\beta)}$ that we use to estimate LAM:
\[LAM \approx 1 - \frac{1}{RR}\sum_{I\in (C_L)} \mathbf{p}^{(3)}_I.\]
dist = distribution(data, th, 3)
rr = rrate(dist)
lam = laminarity(rr, dist)If the distribution has 8 elements, this function will consider line motifs, which makes the process simpler. In this case, we just need the motif with $I = 2$:
\[LAM \approx 1 - \frac{\mathbf{p}^{(3)}_2}{RR}.\]
dist = distribution(data, th, 3; shape = :line)
rr = rrate(dist)
lam = laminarity(rr, dist)Input:
rr: recurrence rate.[probs]: aVector{Float64}returned by the functiondistribution(...).
Output: returns the laminarity as a Float64.
laminarity([x], threshold::Float64; r::Float64)Estimate the laminarity from a data [x]` using a probability distribution and a RR computed from it.
Input:
[x]: input data.[parameter]: set of parameters used to compute the recurrence microstate distribution.r(kwarg): ratio of the total number of microstates to be sampled for the histogram. (defaultr = 0.05)
Output: returns a Tuple{Float64, Float64}.
lam: laminarity asFloat64.rr: recurrence rate asFloat64.
Disorder
Main.RecurrenceMicrostatesAnalysis.disorder — Function
disorder([x], n::Int; kwords...)Compute disorder of a time series [x], based on the work of Flauzino et al. (2025)[10](read here]).
Input:
[x]: input dataset.n: microstates' size, beingn ∈ {2, 3, 4}.
kwargs:
ε: recurrence threshold used as reference to define the range of thresholds where we look for the maximum disorder.ε_min: minimum value of ε of the range.ε_max: maximum value of ε of the range.ε_range_size: number of elements in the range.
Output: return the disorder value as a Float64.
Main.RecurrenceMicrostatesAnalysis.Disorder.get_memory — Function
get_memory(label::Vector{Vector{Int}})Allocate the necessary memory to compute disorder. Return a Vector{Float64}.
Main.RecurrenceMicrostatesAnalysis.Disorder.get_class_probs! — Function
get_class_probs!(memory::Vector{Float64}, probs::Vector{Float64}, label::Vector{Vector{Int}}, class::Int)Get the probabilities associated to a specific motif' class for a given motif size n. memory is a block of memory allocate by the function get_memory(n) that will be modified to store the probabilities from probs associated to the class.
Main.RecurrenceMicrostatesAnalysis.Disorder.get_norm_class_probs! — Function
get_norm_class_probs!(memory::Vector{Float64}, probs::Vector{Float64}, label::Vector{Vector{Int}}, class::Int)Get the normalized probabilities associated to a specific motif' class for a given motif size n. memory is a block of memory allocate by the function get_memory(n) that will be modified to store the probabilities from probs associated to the class.
Utilitary Functions
Main.RecurrenceMicrostatesAnalysis.prepare — Function
prepare([solution], σ::Union{Float64, Int}; transient::Int, K::Int)Prepare a problem solved by the library DifferentialEquations.jl to be used in RecurrenceMicrostatesAnalysis.jl. This function applies the sampling parameter (σ) to discretize the continuous time series, as proposed by [11].
Input:
[solution]: solution returned by the libraryDifferentialEquations.jl.σ: sampling parameter; it defines the time resolution of discretized data.transient(kwarg): number of points, without application of sampling, that will be ignored.K(kwarg): maximum length of the returned data series.
Output:
data: returns the prepared data in the format of aMatrix{Float64}. Each row represents a system component, and each column represents a time step.
Main.RecurrenceMicrostatesAnalysis.find_parameters — Function
find_parameters([x], n::Int; r::Float64 = 0.05, ε_max_range = 0.5)This function calculates the maximum microstate entropy for the RP of the input time series given the microstate size n and the input ratio r of the total number of microstates.
Input:
[x]: input data.n: microstate size.r(kwarg): ratio of the total number of microstates to be sampled for the histogram (defaultr = 0.05)ε_max_range(kwarg): percentage of the maximum distance to be used as the range in the process. (defaultε_max_range = 0.5)fraction(kwarg): iteration fraction. (defaultfraction = 5)shape(kwarg): motif shape. (defaultshape = :square)
Output: (is a Tuple{Float64, Float64})
εopt: the value of the vicinity parameter that maximizes the recurrence microstates entropy.Smax: maximum recurrence microstates entropy for the input time series.
Recurrence Functions
Main.RecurrenceMicrostatesAnalysis.recurrence — Function
recurrence([x], [y], parameters, idx::AbstractVector{Int}, metric::Metric, dim::AbstractVector{Int})Compute the recurrence between two position defined by idx of the datasets x and y. parameters defines which type of recurrence it will use, like the standard recurrence and the recurrence with corridor threshold. metric defines the norm applied to the datasets, it uses the library Distances.jl to compute the metric. dim is just used for high-dimensional problems, such as images. A recurrence function is a function of the form
\[R_{ij} = \Theta(\varepsilon - \|\mathbf{x}_i - \mathbf{y}_j\|),\]
where $\Theta$ is the Heaviside function, $\|\cdot\|$ denotes an appropriate norm, and $\varepsilon$ is a threshold parameter that defines the maximum distance between two points for them to be considered $\varepsilon$-recurrent to each other.
# Examples
# For a 2D recurrence space.
@inline function recurrence(...)
return @inbounds evaluate(metric, x[idx[1]], y[idx[2]]) <= threshold
end
# For a recurrence tensor space (generalization to spatial data)
@inline function recurrence(...)
return @inbounds evaluate(metric, view(x, :, view(idx, 1:dim[1])), view(y, :, view(idx, dim[1]+1:dim[1] + dim[2]))) <= threshold
endInput:
[x]: a dataset.[y]: a dataset.[parameter]: set of parameters used to compute the recurrence microstate distribution, i.e., the value of $\varepsilon$.[idx]: vector of indeces from[x]and[y]to calculate the recurrence between them.metric: metric fromDistances.jlused to compute the recurrence. It defines the norm $\|\cdot\|$.[dim]: number of dimensions derived from[x]and[y]. If you are using a time series it is usually[1,1].
Output: Recurrence functions return true when we have a recurrence, and false otherwise.