Type: | Package |
Title: | Unidimensional and Multidimensional Methods for Recurrence Quantification Analysis |
Version: | 2.0.7 |
Date: | 2025-06-25 |
Maintainer: | Moreno I. Coco <moreno.cocoi@gmail.com> |
Description: | Auto, Cross and Multi-dimensional recurrence quantification analysis. Different methods for computing recurrence, cross vs. multidimensional or profile iti.e., only looking at the diagonal recurrent points, as well as functions for optimization and plotting are proposed. in-depth measures of the whole cross-recurrence plot, Please refer to Coco and others (2021) <doi:10.32614/RJ-2021-062>, Coco and Dale (2014) <doi:10.3389/fpsyg.2014.00510> and Wallot (2018) <doi:10.1080/00273171.2018.1512846> for further details about the method. |
Depends: | R (≥ 3.5.0) |
Imports: | Matrix, pracma, rdist, tseriesChaos, gplots, dplyr, ggplot2 |
License: | GPL (≥ 3) |
Collate: | 'crqa.R' 'crqa_helpers.R' 'drpfromts.R' 'lorenzattractor.R' 'mdDelay.R' 'mdFnn.R' 'optimizeParam.R' 'piecewiseRQA.R' 'plot_rp.R' 'simts.R' 'spdiags.R' 'wincrqa.R' 'windowdrp.R' |
RoxygenNote: | 7.3.2 |
LazyData: | true |
VignetteBuilder: | knitr |
Encoding: | UTF-8 |
Suggests: | knitr, rmarkdown |
NeedsCompilation: | yes |
Packaged: | 2025-06-25 09:46:08 UTC; moreno |
Author: | Moreno I. Coco [cre, aut], Dan Mønster [aut], Giuseppe Leonardi [aut], Rick Dale [aut], Sebastian Wallot [aut], James D. Dixon [ctb], John C. Nash [ctb], Alexandra Paxton [ctb] |
Repository: | CRAN |
Date/Publication: | 2025-06-25 10:10:01 UTC |
Unidimensional and Multidimensional Methods for Recurrence Quantification Analysis
Description
Auto, Cross and Multi-dimensional recurrence quantification analysis. Different methods for computing recurrence, cross vs. multidimensional or profile iti.e., only looking at the diagonal recurrent points, as well as functions for optimization and plotting are proposed. in-depth measures of the whole cross-recurrence plot, Please refer to by Coco and Dale (2014) <doi:10.3389/fpsyg.2014.00510> and Wallot (2018) <doi: 10.1080/00273171.2018.1512846> for further details about the method.
Details
crqa
: Core recurrence function, which examines
recurrent structures of a single rqa, two crqa,
or multidimensional time-series mdcrqa, which are
time-delayed and embedded in higher dimensional space.
The approach compares the phase space trajectories of
the time-series in the same phase-space when delays are
introduced. A distance matrix between the time-series,
delayed and embedded is calculated. Several measures
representative of the underlying dynamics of the system
are extracted.
drpfromts
: Method to explore the diagonal profile of
the recurrence plot (Auto, Cross, or Multi-dimensional).
It returns the recurrence for different delays,
the maximal recurrence observed and the delay at which it occurred.
lorenzattractor
: An implementation of the Lorenz dynamical system,
which describes the motion of a possible particle, which will
neither converge to a steady state, nor diverge to infinity;
but rather stay in a bounded but 'chaotically' defined
region, i.e., an attractor.
mdDelay
:Estimates time delay for embedding of a
multi-dimensional dataset.
mdFnn
: Computes the percentage of false nearest
neighbors for multidimensional time series as a function
of embedding dimension.
optimizeParam
: Iterative procedure to examine the values
of delay, embedding dimension and radius
to compute recurrence plots of one, two,
or more time-series.
piecewiseRQA
: This is a convenience function which breaks down
the computation of large recurrence plots into a collection of
smaller recurrence plots. It can ease speed and memory issues
if an appropriate size for the block is found.
plotRP
: A convenience function to plot the RP matrix returned by
the crqa.
simts
: A simple algorithm for producing a time-series that drives
a second time-series using parameters, which change independent and conditional
probability of an event to occur.
wincrqa
: A recurrence plot is computed in overlapping
windows of a certain size for a number of delays smaller
than the size of the window; and measures of it extracted.
windowdrp
: A recurrence plot is computed in overlapping
windows of a specified size for a number of
delays smaller than the size of the window.
In every window, the recurrence value for the
different delays is calculated.
A mean is then taken across the delays to obtain
a recurrence value in that particular window.
Author(s)
Moreno I. Coco moreno.cocoi@gmail.com Dan Monster danm@econ.au.dk Giuseppe Leonardi g.leonardi@vizja.pl Rick Dale rdale@ucla.edu Sebastian Wallot sebastian.wallot@ae.mpg.de
References
Webber Jr, C. L., and Zbilut, J. P. (2005). Recurrence quantification analysis of nonlinear dynamical systems. Tutorials in contemporary nonlinear methods for the behavioral sciences, 26-94.
Marwan, N., and Kurths, J. Nonlinear analysis of bivariate data with cross recurrence plots. Physics Letters A 302.5 (2002): 299-307.
Coco, M. I., Monster, D., Leonardi, G., Dale, R., & Wallot, S. (2021). Unidimensional and Multidimensional Methods for Recurrence Quantification Analysis with crqa. R Journal, 13(1).
Examples
# use the available data
data(crqa)
listener = eyemovement$listener
narrator = eyemovement$narrator
delay = 1; embed = 1; rescale = 0; radius = .1;
normalize = 0; mindiagline = 2; minvertline = 2;
tw = 0; whiteline = FALSE; recpt = FALSE; side = "both"
method = 'crqa'; metric = 'euclidean';
datatype = "categorical"
ans = crqa(narrator, listener, delay, embed, rescale, radius, normalize,
mindiagline, minvertline, tw, whiteline, recpt, side, method, metric,
datatype)
print(ans[1:10]) ## last argument of list is the recurrence plot
Eye-movement categorical time-series
Description
Three time series (1250 observations) organised as columns and simulating data from: a periodical sine wave, one of the dimensions of the Lorenz attractor and a white noise signal.
Usage
Figure_1
Format
A data frame with 1250 rows and 3 columns:
- sinus
A periodical sine wave
- lorenz
One dimension of a Lorenz attractor
- wnoise
White noise
A unidimensional sinusoidal time series
Description
A unidimensional sinusoidal time series
Usage
Figure_2
Format
A data frame with 45 rows and 1 column:
- V1
A unidimensional sinusoidal time series
Simulated time series of the three dimensions from the Lorenz system
Description
Simulated time series of the three dimensions from the Lorenz system
Usage
Figure_3
Format
A data frame with 2048 rows and 3 columns:
- V1
First dimension of the Lorenz system
- V2
Second dimension of the Lorenz system
- V3
Third dimension of the Lorenz system
Figure_6
Description
Speed (time in seconds) and memory (peak RAM in MB) performance of crqa() and piecewiseRQA() for simulated data of increasing size (from 3000 to 7000 data points), compared in blocks of different sizes (from 1000 to 6500 in increments of 500)
Usage
Figure_6
Format
A data frame with 28 rows and 2 variables:
- speed
The time in seconds to perform a crqa() analysis
- memory
The peak MB occupied in the RAM to perform a crqa() analysis
- typeRQA
Whether the analysis contained all data points (full) or was computed piecewise (piece)
- datapoint
The number of datapoints
- blocksize
The dimension of the block
Auto, cross and multidimensional recurrence measures of one, two or multiple time-series, time-delayed and embedded in higher dimensional space
Description
Core recurrence function, which examines recurrent structures of a single (rqa), two (crqa), or multidimensional time-series (mdcrqa), which are time-delayed and embedded in higher dimensional space. The approach compares the phase space trajectories of the time-series in the same phase-space when delays are introduced. A distance matrix between the time-series, delayed and embedded is calculated. Several measures representative of the underlying dynamics of the system are extracted (explained below).
Usage
crqa(ts1, ts2, delay, embed, rescale, radius, normalize,
mindiagline, minvertline, tw, whiteline, recpt, side, method,
metric, datatype)
Arguments
ts1 |
First time-series dataset. |
ts2 |
Second time-series dataset |
delay |
The delay unit by which the series are lagged. |
embed |
The number of embedding dimension for phase-reconstruction, i.e., the lag intervals. |
rescale |
Rescale the distance matrix; if rescale = 0 (do nothing); if rescale = 1 (mean distance of entire matrix); if rescale = 2 (maximum distance of entire matrix). if rescale = 3 (minimum distance of entire matrix). if rescale = 4 (euclidean distance of entire matrix). |
radius |
A threshold, cut-off, constant used to decide whether two points are recurrent or not. |
normalize |
Normalize the time-series; if normalize = 0 (do nothing); if normalize = 1 (Unit interval); if normalize = 2 (z-score). |
mindiagline |
A minimum diagonal length of recurrent points. Usually set to 2, as it takes a minimum of two points to define any line. |
minvertline |
A minimum vertical length of recurrent points. |
tw |
The Theiler window parameter |
whiteline |
A logical flag to calculate (TRUE) or not (FALSE) empty vertical lines. |
recpt |
A logical flag indicating whether measures of cross-recurrence are calculated directly from a recurrent plot (TRUE) or not (FALSE). |
side |
A string indicating whether recurrence measures should be calculated in the 'upper' triangle of the RP 'lower' triangle of the matrix, or 'both'. LOC is automatically excluded for 'upper' and 'lower'. |
method |
A string to indicate the type of recurrence analysis to peform. There are three options: rqa (autorecurrence); crqa(cross-recurrence); mdcrqa(multidimensional recurrence). Default value is crqa |
metric |
A string to indicate the type of distance metric used, default is euclidean but see help rdist() to list all other possible metrics. |
datatype |
a string (continuous or categorical) to indicate whether the nature of the data type |
Details
We recommend setting whiteline = FALSE, as the current version of the library does not make use of such information to extract recurrence measures.
Value
If a recurrence plot (RP) can be calculated and hence recurrence observed the function will returna a list with different measures extracted. Otherwise, the values for the output arguments will be either 0 or NA.
RR |
The percentage of recurrent points falling within the specified radius (range between 0 and 100) |
DET |
Proportion of recurrent points forming diagonal line structures. |
NRLINE |
The total number of lines in the recurrent plot |
maxL |
The length of the longest diagonal line segment in the plot, excluding the main diagonal |
L |
The average length of line structures |
ENTR |
Shannon information entropy of diagonal line lengths longer than the minimum length |
rENTR |
Entropy measure normalized by the number of lines observed in the plot. Handy to compare across contexts and conditions |
LAM |
Proportion of recurrent points forming vertical line structures |
TT |
The average length of vertical line structures |
catH |
Entropy of categorical recurrence plots based on rectangular block structures |
max_vertlength |
The maximum vertical line length |
RP |
The Recurrence Plot sparse matrix data |
Note
Original bits of this code were translated from a Matlab version provided by Rick Dale, and created during the Non-Linear Methods for Psychological Science summer school held at the University of Cincinnati in 2012. The multi-dimensional method for the crqa function has been written together with Sebastian Wallot (sebastian.wallot at aesthetics.mpg.de )
Author(s)
Moreno I. Coco (moreno.cocoi@gmail.com)
References
Coco, M. I., and Dale, R. (2014). Cross-recurrence quantification analysis of categorical and continuous time series: an R package. Frontiers in psychology, 5, 510.
Wallot, S. (2018). Multidimensional Cross-Recurrence Quantification Analysis (MdCRQA) a method for quantifying correlation between multivariate time-series. Multivariate behavioral research, 1-19
See Also
Examples
# use the available data
data(crqa)
listener = eyemovement$listener
narrator = eyemovement$narrator
delay = 1; embed = 1; rescale = 0; radius = .1;
normalize = 0; mindiagline = 2; minvertline = 2;
tw = 0; whiteline = FALSE; recpt = FALSE; side = "both"
method = 'crqa'; metric = 'euclidean';
datatype = "categorical"
ans = crqa(narrator, listener, delay, embed, rescale, radius, normalize,
mindiagline, minvertline, tw, whiteline, recpt, side, method,
metric, datatype)
print(ans[1:11]) ## last argument of list is the recurrence plot
Diagonal recurrence profile
Description
Method to explore the diagonal profile of the recurrence plot (Auto, Cross, or Multi-dimensional). It returns the recurrence for different delays, the maximal recurrence observed and the delay at which it occurred.
Usage
drpfromts(ts1, ts2, windowsize, radius,
delay, embed, rescale, normalize, mindiagline, minvertline, tw,
whiteline, recpt, side, method, metric, datatype)
Arguments
ts1 |
First time-series |
ts2 |
Second time-series |
windowsize |
A constant indicating the range of delays (positive and negative) to explore |
radius |
A threshold, cut-off, constant used to decide whether two points are recurrent or not. |
delay |
The delay unit by which the series are lagged. |
embed |
The number of embedding dimension for phase-reconstruction, i.e., the lag intervals. |
rescale |
Rescale the distance matrix; if rescale = 0 (do nothing); if rescale = 1 (mean distance of entire matrix); if rescale = 2 (maximum distance of entire matrix). if rescale = 3 (minimum distance of entire matrix). if rescale = 4 (euclidean distance of entire matrix). |
normalize |
Normalize the time-series; if normalize = 0 (do nothing); if normalize = 1 (Unit interval); if normalize = 2 (z-score). |
mindiagline |
A minimum diagonal length of recurrent points. Usually set to 2, as it takes a minimum of two points to define any line. |
minvertline |
A minimum vertical length of recurrent points. |
tw |
The Theiler window parameter |
whiteline |
A logical flag to calculate (TRUE) or not (FALSE) empty vertical lines. |
recpt |
A logical flag indicating whether measures of cross-recurrence are calculated directly from a recurrent plot (TRUE) or not (FALSE). |
side |
A string indicating whether recurrence measures should be calculated in the 'upper' triangle of the RP 'lower' triangle of the matrix, or 'both'. LOC is automatically excluded for 'upper' and 'lower'. |
method |
A string to indicate the type of recurrence analysis to peform. There are three options: rqa (autorecurrence); crqa(cross-recurrence); mdcrqa(multidimensional recurrence). Default value is crqa |
metric |
A string to indicate the type of distance metric used, default is euclidean but see help rdist() to list all other possible metrics. |
datatype |
a string (continuous or categorical) to indicate whether the nature of the data type |
Value
A list with the following arguments:
profile |
A vector of recurrence (ranging from 0,1) with length equal to the number of delays explored |
maxrec |
Maximal recurrence observed between the two-series |
maxlag |
Delay at which maximal recurrence is observed |
Author(s)
Moreno I. Coco (moreno.cocoi@gmail.com)
See Also
Examples
# use the available data
data(crqa)
listener = eyemovement$listener
narrator = eyemovement$narrator
res = drpfromts(narrator, listener, windowsize = 100,
radius = 0.001, delay = 1, embed = 1, rescale = 0,
normalize = 0, mindiagline = 2, minvertline = 2,
tw = 0, whiteline = FALSE, recpt = FALSE,
side = 'both', method = 'crqa',
metric = 'euclidean', datatype = 'continuous')
profile = res$profile
plot(seq(1,length(profile),1), profile, type = "l", lwd = 5,
xaxt = "n", xlab = "Lag", ylab = "Recurrence")
Eye-movement categorical time-series
Description
A two-columns dataset of eye-movement fixation scan-pattern, which are temporal sequences of fixated objects.
Usage
eyemovement
Format
A data frame with 1000 rows and 2 variables:
- listener
the listener time series
- narrator
the narrator time series
References
Richardson, D. C., and Dale, R. (2005). Looking to understand: The coupling between speakers and listeners eye movements and its relationship to discourse comprehension. Cognitive Science, 29, 39-54.
Continuous series of hand movements
Description
Hand-movement velocity profiles of two participants (P1 and P2) for the dominant (d) and non-dominant (n) hand.
Usage
handmovement
Format
A dataframe of 5799 observations.
- P1_TT_d
Participant 1 dominant hand
- P1_TT_n
Participant 1 non-dominant hand
- P2_TT_d
Participant 2 dominant hand
- P2_TT_n
Participant 2 non-dominant
References
Wallot, S., Mitkidis, P., McGraw, J. J. and Roepstorff, A. (2016). Beyond synchrony: joint action in a complex production task reveals beneficial effects of decreased interpersonal synchrony. PloS one, 11(12), e0168306.
Simulate the Lorenz Attractor
Description
An implementation of the Lorenz dynamical system, which describes the motion of a possible particle, which will neither converge to a steady state, nor diverge to infinity; but rather stay in a bounded but 'chaotically' defined region, i.e., an attractor.
Usage
lorenzattractor(numsteps, dt, sigma, r, b)
Arguments
numsteps |
The number of simulated points |
dt |
System parameter |
sigma |
System parameter |
r |
System parameter |
b |
System parameter |
Value
It returns a matrix with the 3 dimensions of the Lorenz
Author(s)
Moreno I. Coco (moreno.cocoi@gmail.com)
References
Lorenz, Edward Norton (1963). Deterministic nonperiodic flow. Journal of the Atmospheric Sciences 20(2) 130-141.
Examples
## initialize the parameters
numsteps = 2 ^ 11; dt = .01; sigma = 10; r = 28; b = 8/3;
res = lorenzattractor(numsteps, dt, sigma, r, b)
Find optimal delay from a multi-dimensional dataset.
Description
Estimates time delay for embedding of a multi-dimensional dataset.
Usage
mdDelay(data, nbins, maxlag, criterion, threshold)
Arguments
data |
The matrix containing all variables |
nbins |
The number of bins considered to estimate mutual information |
maxlag |
Number of lags considered |
criterion |
A string to indicate what delay optizes mutual information: 'firstBelow' uses the lowest delay at which the AMI function drops below the value set by the threshold parameter. 'localMin' uses the position of the first local minimum of the AMI function. The categorical state on which phi is calculated |
threshold |
Value to select the delay when AMI drops below it. |
Value
It returns the recurrence phi-coefficient profile for state k for all delays considered
Author(s)
Sebastian Wallot, Max Planck Insitute for Empirical Aesthetics Dan Moenster, Aarhus University, Moreno I. Coco, University of East London
References
Wallot, S., and Moenster, D. (2018). Calculation of average mutual information (AMI) and false-nearest neighbors (FNN) for the estimation of embedding parameters of multidimensional time -series in Matlab. Front. Psychol. - Quantitative Psychology and Measurement
See Also
Examples
nbins = 10; maxlag = 10; criterion = "firstBelow"; threshold = exp(-1)
data(crqa) ## load the data
handset = handmovement[1:300, ] ## take less points
mdDelay(handset, nbins, maxlag, criterion, threshold)
Find optimal embedding dimension of a multi-dimensional dataset.
Description
Computes the percentage of false nearest neighbors for multidimensional time series as a function of embedding dimension.
Usage
mdFnn(data, tau, maxEmb, numSamples, Rtol, Atol)
Arguments
data |
The matrix of data to estimate FNN. |
tau |
Time delay for embedding. |
maxEmb |
Maximum number of embedding dimensions considered |
numSamples |
Number of randomly drawn coordinates from phase-space used to estimate FNN |
Rtol |
First distance criterion for separating false neighbors |
Atol |
Second distance criterion for separating false neighbors |
Value
It returns the percentage of false neighbors for each embedding.
Author(s)
Sebastian Wallot, Max Planck Insitute for Empirical Aesthetics Dan Moenster, Aarhus University, Moreno I. Coco, University of East London
References
Kennel, M. B., Brown, R., & Abarbanel, H. D. (1992). Determining embedding dimension for phase-space reconstruction using a geometrical construction. Physical review A, 45, 3403. Wallot, S., and Moenster, D. (2018). Calculation of average mutual information (AMI) and false-nearest neighbors (FNN) for the estimation of embedding parameters of multidimensional time-series in Matlab. Front. Psychol. - Quantitative Psychology and Measurement
See Also
Examples
tau = 1; maxEmb = 10; numSamples = 500; Rtol = 10; Atol = 2
data(crqa) ## load the data
handset = handmovement[1:300, ] ## take less points
mdFnn(handset, tau, maxEmb, numSamples, Rtol, Atol)
Estimate optimal delay, embedding dimension and radius for continuous time-series data
Description
Iterative procedure to examine the values of delay, embedding dimension and radius to compute recurrence plots of one, two, or more time-series.
Usage
optimizeParam(ts1, ts2, par, min.rec, max.rec)
Arguments
ts1 |
First time-series |
ts2 |
Second time-series |
par |
A list of parameters needed for the optimization, refer to the Details section. |
min.rec |
The minimum value of recurrence accepted. Default = 2 |
max.rec |
The maximum value of recurrence accepted. Default = 5 |
Details
The optimization can be applied both to uni-dimensional time-series (method = crqa), or multi-dimensional (method = mdcrqa)
The procedure is identical in both cases:
1) Identify a delay that accommodates both time-series by finding the local minimum where mutual information between them drops, and starts to level off. When one ts has a considerably longer delay indicated than the another, the function selects the longer delay of the two to ensure that new information is gained for both. When the delays are close to each other, the function computes the mean of the two delays.
2) Determine embedding dimensions by using false nearest neighbors and checking when it bottoms out (i.e., there is no gain in adding more dimensions). If the embedding dimension for the two ts are different the algorithm selects the higher embedding dimension of the two to make sure that both time series are sufficiently unfolded.
3) Determine radius yielding a recurrence rate between 2-5 To do so, we first determine a starting radius that yields approximately 25 We generate a sampled sequence of equally spaced possible radi from such radius till 0, using as unit for the sequence step, the standard deviation of the distance matrix divided by a scaling parameter (radiusspan). The larger this parameter, the finer the unit.
For uni-dimensional time-series, the user has to decide how to
choose the value of average mutual information (i.e., typeami = mindip
,
the lag at which minimal information is observed, or typeami = maxlag
,
the maximum lag at which minimal information is observed) and the relative
percentage of information gained in FNN, relative to the first embedding
dimension, when higher embeddings are considered (i.e., fnnpercent
).
Then, as crqa
is integrated in the optimizeParam
to estimate the
radius, most of the arguments are the same (e.g., mindiagline
or tw
).
For multidimensional series, the user needs to specify the right RQA method
(i.e., method = "mdcrqa"
). Then, for the estimation of the delay via
AMI: (1) nbins
the number of bins to compute the two-dimensional
histogram of the original and delayed time series and (2) the criterion
to select the delay (firstBelow
to use the lowest delay at which the AMI
function drops below the value set by the threshold
argument, and
localMin
to use the position of the first local AMI minimum).
The estimation of the embedding dimensions instead needs the following arguments:
(1) maxEmb
, which is the maximum number of embedding dimensions considered,
(2) noSamples
, which is the number of randomly drawn coordinates
from phase-space used to estimate the percentage of false-nearest neighbors,
(3) Rtol
, which is the first distance criterion for separating false
neighbors, and (4) Atol
, which is the second distance criterion for
separating false neighbors. The radius is estimated as before.
Value
It returns a list with the following arguments:
radius |
The optimal radius value found |
emddim |
Number of embedding dimensions |
delay |
The lag parameter. |
Note
As optimizeParam
uses crqa
to estimate
the parameters: the additional arguments normalize,
rescale, mindiagline, minvertline, whiteline,
recpt should be supplied in the par list.
Set up relatively large radiusspan (e.g. 100),
for a decent coverage of radius values.
Author(s)
Moreno I. Coco (moreno.cocoi@gmail.com), James A. Dixon (james.dixon@uconn.edu) Sebastian Wallot, Max Planck Insitute for Empirical Aesthetics Dan Moenster, Aarhus University
References
Marwan, N., Carmen Romano, M., Thiel, M., and Kurths, J. (2007). Recurrence plots for the analysis of complex systems. Physics Reports, 438(5), 237-329.
See Also
Examples
data(crqa) ## load the data
handset = handmovement[1:300, ] ## take less points
P1 = cbind(handset$P1_TT_d, handset$P1_TT_n)
P2 = cbind(handset$P2_TT_d, handset$P2_TT_n)
par = list(method = "mdcrqa", metric = "euclidean", maxlag = 20,
radiusspan = 100, radiussample = 40, normalize = 0,
rescale = 4, mindiagline = 10, minvertline = 10, tw = 0,
whiteline = FALSE, recpt = FALSE, side = "both",
datatype = "continuous", fnnpercent = NA,
typeami = NA, nbins = 50, criterion = "firstBelow",
threshold = 1, maxEmb = 20, numSamples = 500,
Rtol = 10, Atol = 2)
results = optimizeParam(P1, P2, par, min.rec = 2, max.rec = 5)
print(unlist(results))
Compute recurrence plots for long time-series data series using a block (piece-wise) method.
Description
This is a convenience function which breaks down the computation of large recurrence plots into a collection of smaller recurrence plots. It can ease speed and memory issues if an appropriate size for the block is found.
Usage
piecewiseRQA(ts1, ts2, blockSize, delay, embed, rescale, radius,
normalize, mindiagline, minvertline, tw, whiteline, recpt, side,
method, metric, datatype, typeRQA, windowsize)
Arguments
ts1 |
First time-series. |
ts2 |
Second time-series. |
blockSize |
The dimension of the time-series subunit in which the-recurrence plot will be computed |
delay |
The delay unit by which the series are lagged. |
embed |
The number of embedding dimension for phase-reconstruction, i.e., the lag intervals. |
rescale |
Rescale the distance matrix; if rescale = 0 (do nothing); if rescale = 1 (mean distance of entire matrix); if rescale = 2 (maximum distance of entire matrix). if rescale = 3 (minimum distance of entire matrix). if rescale = 4 (euclidean distance of entire matrix). |
radius |
A threshold, cut-off, constant used to decide whether two points are recurrent or not. |
normalize |
Normalize the time-series; if normalize = 0 (do nothing); if normalize = 1 (Unit interval); if normalize = 2 (z-score). |
mindiagline |
A minimum diagonal length of recurrent points. Usually set to 2, as it takes a minimum of two points to define any line. |
minvertline |
A minimum vertical length of recurrent points. |
tw |
The Theiler window parameter |
whiteline |
A logical flag to calculate (TRUE) or not (FALSE) empty vertical lines. |
recpt |
A logical flag indicating whether measures of cross-recurrence are calculated directly from a recurrent plot (TRUE) or not (FALSE). |
side |
A string indicating whether recurrence measures should be calculated in the 'upper' triangle of the RP 'lower' triangle of the matrix, or 'both'. LOC is automatically excluded for 'upper' and 'lower'. |
method |
A string to indicate the type of recurrence analysis to peform. There are three options: rqa (autorecurrence); crqa(cross-recurrence); mdcrqa(multidimensional recurrence). Default value is crqa |
metric |
A string to indicate the type of distance metric used, default is euclidean but see help rdist() to list all other possible metrics. |
datatype |
a string (continuous or categorical) to indicate whether the nature of the data type |
typeRQA |
a string (full or diagonal) to indicate whether piecewise recurrence quantification measures should be returned for full plot or for the diagonal profile |
windowsize |
the size of the window around the diagonal of the recurrence (if typeRQA = diagonal) |
Details
It is important to estimate the size of the block that may maximize the speed of the computation. We suggest to monitor how speed and memory usage changes as a result of increasing the time-series and the block size. We also recommend setting whiteline = FALSE, as the current version of the library does not make use of such information to extract measures of recurrence.
Value
If an RP can be calculated and recurrence is found, the piecewiseRQA will return exactly the same measures as crqa() if the typeRQA is set to 'full' and drpdfromts() if the typeRQA is set to 'diagonal'. Please refer to the help file for those two functions for details about the measures.
RP |
The Recurrence Plot sparse matrix data |
Author(s)
Moreno I. Coco (moreno.cocoi@gmail.com) based on Matlab code by Sebastian Wallot
See Also
Examples
## Uncomment and run locally
## generate some data using pracma
# ts1 = seq(0.1, 200, .1)
# ts1 = sin(ts1) + linspace(0, 1,length(ts1));
# ts2 = ts1
## initialize the parameters
# blockSize = 100; delay = 15; embed = 2; rescale = 0; radius = 0.5;
# normalize = 0; mindiagline = 2; minvertline = 2;
# tw = 1; whiteline = FALSE; recpt = FALSE; side = 'both'
# method = "crqa"; metric = 'euclidean'; datatype = "continuous"
# typeRQA = "full"; windowsize = NA
# pieceRP = piecewiseRQA(ts1, ts2, blockSize, delay, embed, rescale,
# radius, normalize, mindiagline, minvertline,
# tw, whiteline, recpt, side,
# method, metric, datatype, typeRQA,
# windowsize)
# print(unlist(pieceRP[1:10]))
Plot a recurrence matrix
Description
Given a recurrence plot object (RP) returned from the function crqa()
,
visualize it.
Usage
plot_rp(rp_matrix,
title = "",
xlabel = "Time",
ylabel = "Time",
pcolour = "black",
geom = c("tile", "point", "void"),
flip_y = FALSE)
Arguments
rp_matrix |
A recurrence plot sparse matrix from |
title |
Main title of the plot (character). Default: |
xlabel |
The x-axis label (character). Default: |
ylabel |
The y-axis label (character). Default: |
pcolour |
The colour used for points in the plot (character). Default: |
geom |
The ggplot2 geom used to draw the points (character).
One of: |
flip_y |
Logical argument controlling whether to flip the directionality of
the y axis (logical). Default: |
Details
The argument geom
is a character denoting what type of ggplot2 geom is
used to plot the points in the recurrence plot. Allowed values are
"tile", "point"
and "void"
. The latter produces an empty plot,
so the user is required to add their own geom to the output. See example below.
Value
A ggplot2 plot.
Author(s)
Mønster, D. danm@econ.au.dk
See Also
This plotting solution as substituted what was previously plotRP()
available till version crqa.2.0.5
Examples
# Here is a very short example based on these two time series
ts_1 <- c(0, 0, 1, 1, 0, 0, 2, 2, 1, 1)
ts_2 <- c(1, 1, 2, 2, 0, 0, 1, 2, 2, 1)
# Create the cross recurrence matrix
short_rec <- crqa(ts_1, ts_2,
delay = 1, embed = 1, radius = 0.001,
rescale = 1, method = "crqa")
# Extract the cross recurrence matrix
small_rp <- short_rec$RP
# Make a cross recurrence plot with blue tiles
plot_rp(small_rp, pcolour = "blue", geom = "tile")
# Make a cross recurrenbce plot with blue points
plot_rp(small_rp, pcolour = "blue", geom = "point")
## Uncomment below if you want to have more example. Not runnable because
## of CRAN, globally unavailable function plus timing constraint to check vignette.
# Neither of the two plots above may be suitable
# Using geom = "void" we can add a custumized geom, here blue tiles
# that are smaller than the default width and height
# plot_rp(small_rp, geom = "void") +
# geom_tile(fill = "blue", height = 0.75, width = 0.75)
# Another custom geom uses tiles with a border colour
# plot_rp(small_rp, geom = "void") +
# geom_tile(fill = "blue", colour = "pink", linewidth = 1)
# Use the eyemovement dataset to create a cross-recurrence plot
# large_crqa <- crqa(eyemovement$narrator, eyemovement$listener,
# delay = 1, embed = 1, radius = 0.001,
# method = "crqa", metric = "euclidean", datatype = "continuous")
# Extract the recurrence matrix from the result
# large_rp <- large_crqa$RP
# Create a recurrence plot using defaults
# plot_rp(large_rp)
# The same recurrence plot with flipped y axis and using geom_point()
# plot_rp(large_rp, flip_y = TRUE, geom = "point", title = "Flipped y axis")
# Add axes labels, a title and extra plot elements
# Note that we can add to the output, here added the line of synchonization
# using ggplot2's geom_abline() function.
# plot_rp(large_rp,
# xlabel = "Narrator",
# ylabel = "Listener",
# title = "Coloured tiles with line of synchronization",
# pcolour = "blue") +
# geom_abline(intercept = 0, slope = 1)
Simulate dichotomous binary time-series
Description
A simple algorithm for producing a time-series that drives a second time-series (1 for event occurrence; 0 otherwise) using parameters, which change independent and conditional probability of an event to occur.
Usage
simts(BL1, BL2, BLR1, BLR2, BL2C1, tsL)
Arguments
BL1 |
Base event rate of the first time-series |
BL2 |
Base event rate of the second time-series |
BLR1 |
Rate of repetition in the first series |
BLR2 |
Rate of repetition in the second series |
BL2C1 |
Conditional probability of repetition. |
tsL |
Length of the simulated time-series |
Value
A matrix with two-rows, where the first row is
the 'driving' time-series and the second row
is the second time-series. The columns are
the number of simulated points as selected
by the argument tsL
.
Author(s)
Rick Dale and Moreno I. Coco (moreno.cocoi@gmail.com)
Examples
## set up parameters
BL1 = .08; BL2 = .05; BLR1 = .5; BLR2 = .5;
BL2C1 = .33; tsL = 100
ts = simts(BL1, BL2, BLR1, BLR2, BL2C1, tsL)
Extract diagonal matrices
Description
Extracts all nonzero diagonals from the m-by-n matrix A. B is a min(m,n)-by-p matrix whose columns are the p nonzero diagonals of A.
Usage
spdiags(A)
Arguments
A |
An m-by-n matrix with nonzero elements located on p diagonals. |
Details
Compared to the original Matlab implementation: 1) it does not handle the case with more than one input, and 2) (m > n) matrices give the B matrix columns in a different order, but the d vector of indices will also be changed accordingly, so the set of columns is OK, just ordered differently
Value
B |
A min(m,n)-by-p matrix, usually (but not necessarily) full, whose columns are the diagonals of A. |
d |
A vector of length p whose integer components specify the diagonals in A. |
Note
For computational efficiency spdiags is actually
computed using a Fortan implementation (jspd.f
)
Author(s)
John C. Nash (nashjc@uottawa.ca)
Examples
dta <- c(0, 5, 0, 10, 0, 0, 0, 0, 6, 0, 11, 0, 3, 0, 0,
7, 0, 12, 1, 4, 0, 0, 8, 0, 0, 2, 5, 0, 0, 9)
A1 <- matrix(dta, nrow=5, ncol=6, byrow=TRUE)
print(A1)
res1 <- spdiags(A1)
print(res1)
Categorical sequence of words
Description
A simple text of the nursery rhyme, \'the wheels on the bus\'
Usage
text
Format
A vector of 120 words
References
Verna Hills (1939). The Wheels on the bus. American Childhood (25), 56
Windowed Recurrence Measures
Description
A recurrence plot is computed in overlapping windows of a certain size for a number of delays smaller than the size of the window; and measures of it extracted.
Usage
wincrqa(ts1, ts2, windowstep, windowsize, delay, embed,
radius, rescale, normalize, mindiagline, minvertline, tw, whiteline,
recpt, side, method, metric, datatype, trend)
Arguments
ts1 |
First time-series |
ts2 |
Second time-series |
windowstep |
Interval by which the window is moved. |
windowsize |
The size of the window |
delay |
The delay unit by which the series are lagged. |
embed |
The number of embedding dimension for phase-reconstruction, i.e., the lag intervals. |
radius |
A threshold, cut-off, constant used to decide whether two points are recurrent or not. |
rescale |
Rescale the distance matrix; if rescale = 0 (do nothing); if rescale = 1 (mean distance of entire matrix); if rescale = 2 (maximum distance of entire matrix). if rescale = 3 (minimum distance of entire matrix). if rescale = 4 (euclidean distance of entire matrix). |
normalize |
Normalize the time-series; if normalize = 0 (do nothing); if normalize = 1 (Unit interval); if normalize = 2 (z-score). |
mindiagline |
A minimum diagonal length of recurrent points. Usually set to 2, as it takes a minimum of two points to define any line. |
minvertline |
A minimum vertical length of recurrent points. |
tw |
The Theiler window parameter |
whiteline |
A logical flag to calculate (TRUE) or not (FALSE) empty vertical lines. |
recpt |
A logical flag indicating whether measures of cross-recurrence are calculated directly from a recurrent plot (TRUE) or not (FALSE). |
side |
A string indicating whether recurrence measures should be calculated in the 'upper' triangle of the RP 'lower' triangle of the matrix, or 'both'. LOC is automatically excluded for 'upper' and 'lower'. |
method |
A string to indicate the type of recurrence analysis to peform. There are three options: rqa (autorecurrence); crqa(cross-recurrence); mdcrqa(multidimensional recurrence). Default value is crqa |
metric |
A string to indicate the type of distance metric used, default is euclidean but see help rdist() to list all other possible metrics. |
datatype |
a string (continuous or categorical) to indicate whether the nature of the data type |
trend |
a boolean (TRUE or FALSE) to indicate whether the TREND should be computed of the system |
Value
It returns a matrix where the rows are the
different windows explored, and the columns
are the recurrence measures observed
in that particular window.
Refer to crqa
for the values returned.
Note
If no-recurrence is found in a window, that window will not be saved, and a message about it will be warned. TREND is implemented following a solution proposed by Norbert Marwan, and translated here in R, for those who have asked him. He, however, warns that this measure might strongly depend on the chosen settings to calculate crq. Relying blindly on such measure can, therefore, produce misleading results. Also, we enabled the possibility to input directly a RP with recpt = T, and so extract the measures. This implies that it will not be the same as by conducting phase-space reconstruction on subsets of the time series. So, please, make sure why you are doing it and why.
Author(s)
Moreno I. Coco (moreno.cocoi@gmail.com) Alexandra Paxton alexandra.paxton@uconn.edu
See Also
Examples
data(crqa)
listener = eyemovement$listener[1:200]
narrator = eyemovement$narrator[1:200]
# NB, the parameters for windowsize and windowstep are large to allow
# faster running time, please set them carefully in your analysis.
delay = 1; embed = 1; rescale = 0; radius = 0.001;
normalize = 0; mindiagline = 2; minvertline = 2;
tw = 0; whiteline = FALSE; recpt = FALSE; side = "both"
method = 'crqa'; metric = 'euclidean';
datatype = "continuous";
windowsize = 100; windowstep = 50
trend = FALSE
ans = wincrqa(listener, narrator, windowstep, windowsize, delay, embed,
radius, rescale, normalize, mindiagline, minvertline,
tw, whiteline, recpt, side, method, metric,
datatype, trend)
## other recurrence measures are available in ans
profile = as.numeric(ans$RR)
plot(profile, type = 'l')
Windowed Recurrence Profile
Description
A recurrence plot is computed in overlapping windows of a specified size for a number of delays smaller than the size of the window. In every window, the recurrence value for the different delays is calculated. A mean is then taken across the delays to obtain a recurrence value in that particular window.
Usage
windowdrp(ts1, ts2, windowstep, windowsize, lagwidth,
radius, delay, embed, rescale, normalize, mindiagline, minvertline,
tw, whiteline, recpt, side, method, metric, datatype)
Arguments
ts1 |
First time-series |
ts2 |
Second time-series |
windowstep |
Interval by which the window is moved. |
windowsize |
The size of the window |
lagwidth |
The number of delays to be considered within the window |
radius |
For numeric time-series, the cutoff distance to accept or reject two-points as recurrent |
delay |
The delay unit by which the series are lagged. |
embed |
The number of embedding dimension for phase-reconstruction, i.e., the lag intervals. |
rescale |
Rescale the distance matrix; if rescale = 0 (do nothing); if rescale = 1 (mean distance of entire matrix); if rescale = 2 (maximum distance of entire matrix). if rescale = 3 (minimum distance of entire matrix). if rescale = 4 (euclidean distance of entire matrix). |
normalize |
Normalize the time-series; if normalize = 0 (do nothing); if normalize = 1 (Unit interval); if normalize = 2 (z-score). |
mindiagline |
A minimum diagonal length of recurrent points. Usually set to 2, as it takes a minimum of two points to define any line. |
minvertline |
A minimum vertical length of recurrent points. |
tw |
The Theiler window parameter |
whiteline |
A logical flag to calculate (TRUE) or not (FALSE) empty vertical lines. |
recpt |
A logical flag indicating whether measures of cross-recurrence are calculated directly from a recurrent plot (TRUE) or not (FALSE). |
side |
A string indicating whether recurrence measures should be calculated in the 'upper' triangle of the RP 'lower' triangle of the matrix, or 'both'. LOC is automatically excluded for 'upper' and 'lower'. |
method |
A string to indicate the type of recurrence analysis to peform. There are three options: rqa (autorecurrence); crqa(cross-recurrence); mdcrqa(multidimensional recurrence). Default value is crqa |
metric |
A string to indicate the type of distance metric used, default is euclidean but see help rdist() to list all other possible metrics. |
datatype |
a string (continuous or categorical) to indicate whether the nature of the data type |
Value
It returns a list of arguments where:
profile |
Time-course windowed recurrence profile |
maxrec |
Maximal recurrence observed along the time-course |
maxlag |
The point where maximal recurrence is observed |
Author(s)
Moreno I. Coco (moreno.cocoi@gmail.com) and Rick Dale (rdale@ucmerced.edu)
References
Boker, S. M., Rotondo, J. L., Xu, M., and King, K. (2002). Windowed cross-correlation and peak picking for the analysis of variability in the association between behavioral time series. Psychological Methods, 7(3), 338.
See Also
Examples
# use the available data
data(crqa)
listener = eyemovement$listener
narrator = eyemovement$narrator
# NB, the parameters for windowsize and windowstep are large to allow
# faster running time, please set them carefully in your analysis.
delay = 1; embed = 1; rescale = 1; radius = 0.001;
normalize = 0; mindiagline = 2; minvertline = 2;
tw = 0; whiteline = FALSE; recpt = FALSE; side = "both"
method = 'crqa'; metric = 'euclidean';
datatype = "continuous"; windowsize = 100;
lagwidth = 10; windowstep = 200
ans = windowdrp(narrator, listener, windowstep, windowsize, lagwidth,
radius, delay, embed, rescale, normalize,
mindiagline, minvertline, tw,
whiteline, recpt, side, method, metric,
datatype)
profile = ans$profile; maxrec = ans$maxrec; maxlag = ans$maxlag
plot(profile, type = 'l')