Recurrence Quantification Analysis
The recurrence microstate analysis allows us to estimate values of typical RQA measures, such as determinism and laminarity, with a good precision, and defines some novel quantifiers. We will demonstate in this page how compute these quantifiers using a uniform distribution as input.
julia> using RecurrenceMicrostatesAnalysis, Distributions # Generate our data =Djulia> data = rand(Uniform(0, 1), 3000);julia> th, s = find_parameters(data, 3)(0.2712540022468859, 5.684788996949903)julia> dist = distribution(data, th, 3);
Recurrence Entropy
To compute the recurrence entropy, it is possible to use the rentropy function that receives a recurrence motif probability distribution.
julia> entr = rentropy(dist)5.681533363657713
The rentropy function has also a keyword argument that can be used to ignore some motifs.
julia> entr = rentropy(dist; ignore_motifs = [1, 512])5.549076236101606
Note that the kword ignore_motifs uses as index the notation of Julia, beginning in 1, instead 0. So, the motif 0 is identified by the number 1.
Recurrence Rate
The recurrence rate (RR) can be computed using a similar method to the recurrence entropy.
julia> rr = rrate(dist)0.4625747910794636
Since the recurrence rate is an estimated measure, it has a small error, how you can check in the following graphic, that displays the relative error between the RR computed by RecurrenceMicrostatesAnalysis.jl and the standard approach. 
Determinism
The determinism (DET) can be computed using a recurrence motifs probability distribution and the recurrence rate. It is important to note that it can be done using two motif constrained shapes: :square or :diagonal.
julia> det = determinism(rr, dist)0.7106182486784647julia> det = determinism(rr, distribution(data, th, 3; shape = :diagonal))0.7120854334475366
Similar to RR, the determinism is a quantifier estimated using recurrence microstates analysis, so it has a small error that is demonstrated in the following figure. 
(i) is the uniform distribution, (ii) is the Lorenz system, (iii) is the Logistic map, (iv) is the Rössler system, and (v) is the Bernoulli shifted generalized.
We implement an way to do it without the need to compute the recurrence distributions.
julia> det = determinism(data, th)(0.7120381138764347, 0.46239059486723033)
When we estimate DET directly using this function overload, the library will automatically use a diagonal motif constrained shape.
Laminarity
The laminarity (LAM) can be computed with a method similar to determinism (DET). It is important to note that it can be done using two motif constrained shapes: :square or :line.
julia> lam = laminarity(rr, dist)0.7300716755838004julia> lam = laminarity(rr, distribution(data, th, 3; shape = :line))0.7310626266081572
In the same way, laminarity has a small error associated to it estimation. You can check it in the next figure. 
(i) is the uniform distribution, (ii) is the Lorenz system, (iii) is the Logistic map, (iv) is the Rössler system, and (v) is the Bernoulli shifted generalized.
We implement an way to do it without the need to compute the recurrence distributions.
julia> lam = laminarity(data, th)(0.7295071316160557, 0.4623579587598279)
When we estimate LAM directly using this function overload, the library will automatically use a line motif constrained shape.
Disorder
Disorder is a very powerful quantifier introduced by Flauzino et al.[10] in 2025. The computation of this quantifier is very simple using RecurrenceMicrostatesAnalysis.jl:
julia> microstate_size = 33julia> Ξ = disorder(data, microstate_size)0.9996368268697874
Read more about disorder here.
The function disorder tries to maximize the disorder for a range of $\varepsilon$. To compute it faster, we recommend defining a base value of $\varepsilon$ and a range around it. For example:
julia> Ξ = disorder(data, microstate_size; ε_min = 0.9 * th, ε_max = 1.1 * th, ε_range_size = 6)0.9996422146561869