AccRQA functions
C/C++
-
void accrqa_print_error(Accrqa_Error *error)
Prints error messages associated with error status.
- Parameters:
error – Error status to print.
-
inline void accrqa_version()
Prints accrqa version.
-
void accrqa_LAM(float *output, float *input_data, size_t data_size, int *tau_values, int nTaus, int *emb_values, int nEmbs, int *vmin_values, int nVmins, float *threshold_values, int nThresholds, Accrqa_Distance distance_type, int calc_ENTR, Accrqa_CompPlatform comp_platform, Accrqa_Error *error)
Calculates LAM, TT, TTmax, ENTR and RR RQA measures from supplied time-series. (float)
Array dimensions are as follows, from slowest to fastest varying:
outputis 3D data cube containing RR values, with shape:[
nTaus,nEmbs,nLmins,nThresholds, 5 ].
input_datais 1D and real-valued, with shape:[
data_size]
tau_valuesis 1D and integer-valued, with shape:[
nTaus]
emb_valuesis 1D and integer-valued, with shape:[
nEmbs]
lmin_valuesis 1D and integer-valued, with shape:[
nEmbs]
threshold_valuesis 1D and real-valued, with shape:[
nThresholds]
- Parameters:
output – Multi-dimensional data cube containing RR values.
input_data – Real-valued array of input time-series samples.
data_size – Number of samples (float) of the time-series.
tau_values – Integer array of delay values.
nTaus – Number of delays.
emb_values – Integer array of embedding values.
nEmbs – Number of embeddings.
vmin_values – Integer array of minimal lengths values.
nVmins – Number of minimal lengths.
threshold_values – Real-valued array (float) of threshold values.
nThresholds – Number of threshold values.
distance_type – Norm used in calculation of distance to the line of identity.
calc_ENTR – Turns calculation of ENTR on (1) and off (0).
comp_platform – Compute platform to use.
error – Error status.
-
void accrqa_LAM(double *output, double *input_data, size_t data_size, int *tau_values, int nTaus, int *emb_values, int nEmbs, int *vmin_values, int nVmins, double *threshold_values, int nThresholds, Accrqa_Distance distance_type, int calc_ENTR, Accrqa_CompPlatform comp_platform, Accrqa_Error *error)
Calculates LAM, TT, TTmax, ENTR and RR RQA measures from supplied time-series. (double)
Array dimensions are as follows, from slowest to fastest varying:
outputis 3D data cube containing RR values, with shape:[
nTaus,nEmbs,nLmins,nThresholds, 5 ].
input_datais 1D and real-valued, with shape:[
data_size]
tau_valuesis 1D and integer-valued, with shape:[
nTaus]
emb_valuesis 1D and integer-valued, with shape:[
nEmbs]
lmin_valuesis 1D and integer-valued, with shape:[
nEmbs]
threshold_valuesis 1D and real-valued, with shape:[
nThresholds]
- Parameters:
output – Multi-dimensional data cube containing RR values.
input_data – Real-valued array of input time-series samples.
data_size – Number of samples (double) of the time-series.
tau_values – Integer array of delay values.
nTaus – Number of delays.
emb_values – Integer array of embedding values.
nEmbs – Number of embeddings.
vmin_values – Integer array of minimal lengths values.
nVmins – Number of minimal lengths.
threshold_values – Real-valued array (double) of threshold values.
nThresholds – Number of threshold values.
distance_type – Norm used in calculation of distance to the line of identity.
calc_ENTR – Turns calculation of ENTR on (1) and off (0).
comp_platform – Compute platform to use.
error – Error status.
-
int accrqa_LAM_output_size_in_elements(int nTaus, int nEmbs, int nVmins, int nThresholds)
Calculates size of LAM output array in number of elements.
- Parameters:
nTaus – Number of delays.
nEmbs – Number of embeddings.
nVmins – Number of minimal lengths.
nThresholds – Number of threshold values.
-
void accrqa_DET(float *output, float *input_data, size_t data_size, int *tau_values, int nTaus, int *emb_values, int nEmbs, int *lmin_values, int nLmins, float *threshold_values, int nThresholds, Accrqa_Distance distance_type, int calc_ENTR, Accrqa_CompPlatform comp_platform, Accrqa_Error *error)
Calculates DET, L, Lmax, ENTR and RR RQA measures from supplied time-series. (float)
Array dimensions are as follows, from slowest to fastest varying:
outputis 3D data cube containing RR values, with shape:[
nTaus,nEmbs,nLmins,nThresholds, 5 ].
input_datais 1D and real-valued, with shape:[
data_size]
tau_valuesis 1D and integer-valued, with shape:[
nTaus]
emb_valuesis 1D and integer-valued, with shape:[
nEmbs]
lmin_valuesis 1D and integer-valued, with shape:[
nEmbs]
threshold_valuesis 1D and real-valued, with shape:[
nThresholds]
- Parameters:
output – Multi-dimensional data cube containing RR values.
input_data – Real-valued array of input time-series samples.
data_size – Number of samples (float) of the time-series.
tau_values – Integer array of delay values.
nTaus – Number of delays.
emb_values – Integer array of embedding values.
nEmbs – Number of embeddings.
lmin_values – Integer array of minimal lengths values.
nLmins – Number of minimal lengths.
threshold_values – Real-valued array (float) of threshold values.
nThresholds – Number of threshold values.
distance_type – Norm used in calculation of distance to the line of identity.
calc_ENTR – Turns calculation of ENTR on (1) and off (0).
comp_platform – Compute platform to use.
error – Error status.
-
void accrqa_DET(double *output, double *input_data, size_t data_size, int *tau_values, int nTaus, int *emb_values, int nEmbs, int *lmin_values, int nLmins, double *threshold_values, int nThresholds, Accrqa_Distance distance_type, int calc_ENTR, Accrqa_CompPlatform comp_platform, Accrqa_Error *error)
Calculates DET, L, Lmax, ENTR and RR RQA measures from supplied time-series. (double)
Array dimensions are as follows, from slowest to fastest varying:
outputis 3D data cube containing RR values, with shape:[
nTaus,nEmbs,nLmins,nThresholds, 5 ].
input_datais 1D and real-valued, with shape:[
data_size]
tau_valuesis 1D and integer-valued, with shape:[
nTaus]
emb_valuesis 1D and integer-valued, with shape:[
nEmbs]
lmin_valuesis 1D and integer-valued, with shape:[
nEmbs]
threshold_valuesis 1D and real-valued, with shape:[
nThresholds]
- Parameters:
output – Multi-dimensional data cube containing RR values.
input_data – Real-valued array of input time-series samples.
data_size – Number of samples (double) of the time-series.
tau_values – Integer array of delay values.
nTaus – Number of delays.
emb_values – Integer array of embedding values.
nEmbs – Number of embeddings.
lmin_values – Integer array of minimal lengths values.
nLmins – Number of minimal lengths.
threshold_values – Real-valued array (double) of threshold values.
nThresholds – Number of threshold values.
distance_type – Norm used in calculation of distance to the line of identity.
calc_ENTR – Turns calculation of ENTR on (1) and off (0).
comp_platform – Compute platform to use.
error – Error status.
-
int accrqa_DET_output_size_in_elements(int nTaus, int nEmbs, int nLmins, int nThresholds)
Calculates size of DET output array in number of elements.
- Parameters:
nTaus – Number of delays.
nEmbs – Number of embeddings.
nLmins – Number of minimal lengths.
nThresholds – Number of threshold values.
-
void accrqa_RR(float *output, float *input_data, size_t data_size, int *tau_values, int nTaus, int *emb_values, int nEmbs, float *threshold_values, int nThresholds, Accrqa_Distance distance_type, Accrqa_CompPlatform comp_platform, Accrqa_Error *error)
Calculates RR RQA measure from supplied time-series. (float)
Array dimensions are as follows, from slowest to fastest varying:
outputis 3D data cube containing RR values, with shape:[
nTaus,nEmbs,nThresholds].
input_datais 1D and real-valued, with shape:[
data_size]
tau_valuesis 1D and integer-valued, with shape:[
nTaus]
emb_valuesis 1D and integer-valued, with shape:[
nEmbs]
threshold_valuesis 1D and real-valued, with shape:[
nThresholds]
- Parameters:
output – Multi-dimensional data cube containing RR values.
input_data – Real-valued array of input time-series samples.
data_size – Number of samples (float) of the time-series.
tau_values – Integer array of delay values.
nTaus – Number of delays.
emb_values – Integer array of embedding values.
nEmbs – Number of embeddings.
threshold_values – Real-valued array (float) of threshold values.
nThresholds – Number of threshold values.
distance_type – Norm used in calculation of distance to the line of identity.
comp_platform – Compute platform to use.
error – Error status.
-
void accrqa_RR(double *output, double *input_data, size_t data_size, int *tau_values, int nTaus, int *emb_values, int nEmbs, double *threshold_values, int nThresholds, Accrqa_Distance distance_type, Accrqa_CompPlatform comp_platform, Accrqa_Error *error)
Calculates RR RQA measure from supplied time-series. (double)
Array dimensions are as follows, from slowest to fastest varying:
outputis 3D data cube containing RR values, with shape:[
nTaus,nEmbs,nThresholds].
input_datais 1D and real-valued, with shape:[
data_size]
tau_valuesis 1D and integer-valued, with shape:[
nTaus]
emb_valuesis 1D and integer-valued, with shape:[
nEmbs]
threshold_valuesis 1D and real-valued, with shape:[
nThresholds]
- Parameters:
output – Multi-dimensional data cube containing RR values.
input_data – Real-valued array of input time-series samples.
data_size – Number of samples (double) of the time-series.
tau_values – Integer array of delay values.
nTaus – Number of delays.
emb_values – Integer array of embedding values.
nEmbs – Number of embeddings.
threshold_values – Real-valued array (double) of threshold values.
nThresholds – Number of threshold values.
distance_type – Norm used in calculation of distance to the line of identity.
comp_platform – Compute platform to use.
error – Error status.
-
int accrqa_RR_output_size_in_elements(int nTaus, int nEmbs, int nThresholds)
Calculates size of RR output array in number of elements.
- Parameters:
nTaus – Number of delays.
nEmbs – Number of embeddings.
nThresholds – Number of threshold values.
-
void accrqa_RP(char *output, float *input_data, size_t data_size, int tau, int emb, float threshold, Accrqa_Distance distance_type, Accrqa_Error *error)
Calculates recurrence plot from supplied time-series. (float)
Array dimensions are as follows, from slowest to fastest varying:
outputis 2D data containing recurrence plot in row-major format. Size of the RP output is N = data_size - (emb - 1)*tau, with shape:[
N,N].
input_datais 1D and real-valued, with shape:[
data_size]
- Parameters:
output – Recurrence plot in row-major format.
input_data – Real-valued array of input time-series samples.
data_size – Number of samples (float) of the time-series.
tau – Integer value of delay.
emb – Integer value of embedding.
threshold – Floating point value of threshold.
distance_type – Norm used in calculation of distance to the line of identity.
error – Error status.
-
void accrqa_RP(char *output, double *input_data, size_t data_size, int tau, int emb, double threshold, Accrqa_Distance distance_type, Accrqa_Error *error)
Calculates recurrence plot from supplied time-series. (double)
Array dimensions are as follows, from slowest to fastest varying:
outputis 2D data containing recurrence plot in row-major format. Size of the RP output isN= data_size - (emb - 1)*tau, with shape:[
N,N].
input_datais 1D and real-valued, with shape:[
data_size]
- Parameters:
output – Recurrence plot in row-major format.
input_data – Real-valued array of input time-series samples.
data_size – Number of samples (double) of the time-series.
tau – Integer value of delay.
emb – Integer value of embedding.
threshold – Floating point value of threshold.
distance_type – Norm used in calculation of distance to the line of identity.
error – Error status.
-
void accrqa_RP(int *output, double *input_data, size_t data_size, int tau, int emb, double threshold, Accrqa_Distance distance_type, Accrqa_Error *error)
-
int64_t accrqa_RP_output_size_in_elements(int64_t input_size, int64_t tau, int64_t emb)
Calculates size of RR output array in number of elements.
- Parameters:
nSamples – Number of input samples.
tau – Delay value.
emb – Embedding value.
Python
- accrqa.RP(input_data, tau, emb, threshold, distance_type)
Calculates recurrence plot from supplied time-series. https://en.wikipedia.org/wiki/Recurrence_quantification_analysis
- Parameters:
input_data – The input time-series.
tau – Integer value of delay.
emb – Integer value of embedding.
threshold – Floating point value of threshold.
distance_type – Norm used to calculate distance. Must be instance of
accrqaDistance().
- Returns:
A numpy NDArray containing of RP values.
- Raises:
TypeError – If delay or embedding is negative or zero.
TypeError – If input_data is not numpy.ndarray.
TypeError – If wrong type of the norm to the line of identity is selected.
RuntimeError – If AccRQA library encounters a problem.
- accrqa.RR(input_data, tau_values, emb_values, threshold_values, distance_type, comp_platform=<accrqa.accrqa_compplatform.accrqaCompPlatform object>, tidy_data=True)
Calculates RR measure from supplied time-series. https://en.wikipedia.org/wiki/Recurrence_quantification_analysis
- Parameters:
input_data – The input time-series.
tau_values – Array of delays.
emb_values – Array of embedding values.
threshold_values – Array of threshold values.
distance_type – Norm used to calculate distance. Must be instance of
accrqaDistance().comp_platform – [Optional] Computational platform to be used. Default is cpu. Must be instance of
accrqaCompPlatform().tidy_data – [Optional] Output data in tidy data format. Requires pandas.
- Returns:
A numpy NDArray containing of RR values with dimensions [number of delays, number of embeddings, number of thresholds].
- Raises:
TypeError – If number of delays, embedding or thresholds is zero length.
TypeError – If input_data is not numpy.ndarray.
TypeError – If wrong type of the distance to the line of identity is selected.
TypeError – If wrong computational platform is selected.
RuntimeError – If AccRQA library encounters a problem.
- accrqa.LAM(input_data, tau_values, emb_values, vmin_values, threshold_values, distance_type, calculate_ENTR=True, comp_platform=<accrqa.accrqa_compplatform.accrqaCompPlatform object>, tidy_data=True)
Calculates DET, L, Lmax, ENTR and RR measures from supplied time-series. https://en.wikipedia.org/wiki/Recurrence_quantification_analysis
- Parameters:
input_data – The input time-series.
tau_values – Array of delays.
emb_values – Array of embedding values.
vmin_values – Array of minimal lengths.
threshold_values – Array of threshold values.
distance_type – Norm used to calculate distance. Must be instance of
accrqaDistance().calculate_ENTR – [Optional] Enable calculation of Vmax and ENTR. Default True.
comp_platform – [Optional] Computational platform to be used. Default is cpu. Must be instance of
accrqaCompPlatform().tidy_data – [Optional] Output data in tidy data format. Requires pandas.
- Returns:
A numpy NDArray containing of vertical line based RQA values with dimensions [number of delays, number of embeddings, number of thresholds, RQA_ID]: RQA ID RQA measure: 0 LAM 1 TT 2 TTmax 3 ENTR 4 RR
- Raises:
TypeError – If number of delays, embedding, minimal lengths or thresholds is zero length.
TypeError – If input_data is not numpy.ndarray.
TypeError – If wrong type of the distance to the line of identity is selected.
TypeError – If wrong computational platform is selected.
RuntimeError – If AccRQA library encounters a problem.
- accrqa.DET(input_data, tau_values, emb_values, lmin_values, threshold_values, distance_type, calculate_ENTR=True, comp_platform=<accrqa.accrqa_compplatform.accrqaCompPlatform object>, tidy_data=True)
Calculates DET, L, Lmax, ENTR and RR measures from supplied time-series. https://en.wikipedia.org/wiki/Recurrence_quantification_analysis
- Parameters:
input_data – The input time-series.
tau_values – Array of delays.
emb_values – Array of embedding values.
lmin_values – Array of minimal lengths.
threshold_values – Array of threshold values.
distance_type – Norm used to calculate distance. Must be instance of
accrqaDistance().calculate_ENTR – [Optional] Enable calculation of Lmax and ENTR. Default True.
comp_platform – [Optional] Computational platform to be used. Default is cpu. Must be instance of
accrqaCompPlatform().tidy_data – [Optional] Output data in tidy data format. Requires pandas.
- Returns:
A numpy NDArray containing of diagonal line based RQA values with dimensions [number of delays, number of embeddings, number of thresholds, RQA_ID]: RQA ID RQA measure: 0 DET 1 L 2 Lmax 3 ENTR 4 RR
- Raises:
TypeError – If number of delays, embedding, minimal lengths or thresholds is zero length.
TypeError – If input_data is not numpy.ndarray.
TypeError – If wrong type of the distance to the line of identity is selected.
TypeError – If wrong computational platform is selected.
RuntimeError – If AccRQA library encounters a problem.
- accrqa.RR_target(input_data, tau, emb, target_RR, distance_type, epsilon=0.01, comp_platform=<accrqa.accrqa_compplatform.accrqaCompPlatform object>, max_iter=20, threshold_min=0, threshold_max=10)
Finds the recurrence rate threshold associated with a target recurrence rate (RR) value using a bisection search algorithm within specified precision.
This function performs a binary search between threshold_min and threshold_max to find the threshold value that yields a recurrence rate closest to the target_RR within the specified epsilon tolerance or maximum number of iterations.
Parameters:
- input_dataNDArray
Input time series data as a numpy array
- tauint
Time delay parameter for phase space reconstruction (must be > 0)
- embint
Embedding dimension for phase space reconstruction (must be > 0)
- target_RRfloat
Target recurrence rate value between 0 and 1
- distance_typeaccrqaDistance
Distance metric instance for recurrence plot calculations
- epsilonfloat, optional
Precision tolerance for convergence (default: 0.01, must be > 0)
- comp_platformaccrqaCompPlatform, optional
Computational platform instance (default: NVIDIA GPU)
- max_iterint, optional
Maximum number of bisection iterations (default: 20, must be > 0)
- threshold_minfloat, optional
Lower bound for threshold search range (default: 0)
- threshold_maxfloat, optional
Upper bound for threshold search range (default: 10)
Returns:
- float
The threshold value that produces a recurrence rate closest to target_RR within the specified epsilon tolerance
- float
The RR for given threshold.
Raises:
- TypeError
If distance_type is not an accrqaDistance instance
If comp_platform is not an accrqaCompPlatform instance
If input_data is not a numpy array
If tau or emb are <= 0
If target_RR is not between 0 and 1
If epsilon <= 0
If max_iter <= 0
If threshold_min >= threshold_max
If threshold_min yields RR > target_RR
If threshold_max yields RR < target_RR
If threshold_min and threshold_max yield identical RR values
- accrqa.accrqaDistance(distance_type)
A class representing a distance type used in AccRQA.
This class maps supported distance type names (e.g., “euclidean”, “maximal”) to integer identifiers.
- Parameters:
distance_type (str) – The name of the distance type to use. Must be one of supported distance types.
- Raises:
ValueError – If the provided distance_type is not supported.
- Supported distance types:
euclidian: L-2 norm or an euclidean norm.
maximal: L-inf or a maximal norm.
- accrqa.accrqaCompPlatform(comp_platform)
A class representing a coputational platform supported by AccRQA.
This class maps supported computational platforms (e.g., “cpu”, “nv_gpu”) to integer identifiers.
- Parameters:
comp_platform (str) – The name of the computational platform to use. Must be one of supported computational platforms.
- Raises:
ValueError – If the provided comp_platform is not supported.
- Supported computational platforms:
cpu: CPU parallelised using OpenMP.
nv_gpu: NVIDIA GPU parallelised using CUDA.
R
- accrqa_DET(input, tau_values, emb_values, lmin_values, threshold_values, norm = "euclidean", calc_ENTR = TRUE, platform)
This function calculates the determinism (DET) for cross-recurrence quantification analysis (CRQA) based on a set of input parameters, including time delay, embedding dimensions, minimum line length, threshold values, and normalization.
The function performs cross-recurrence quantification analysis (CRQA) using the specified parameters.
DET measures the proportion of recurrent points forming diagonal lines in the recurrence plot,
which indicates deterministic structures in the data. If calc_ENTR is TRUE, the entropy of diagonal
line lengths is also computed.
- Parameters:
- inputnumeric matrix or data frame
A numeric matrix or data frame representing the input data for CRQA analysis.
- tau_valuesnumeric vector
A numeric vector specifying the time delay(s) to be used in the analysis.
- emb_valuesnumeric vector
A numeric vector specifying the embedding dimensions to be tested.
- lmin_valuesnumeric vector
A numeric vector specifying the minimum diagonal line lengths for DET computation.
- threshold_valuesnumeric vector
A numeric vector specifying the threshold values for recurrence computation.
- normcharacter string, optional
A character string specifying the normalization method to be used. Options may include
"euclidean","maximal", etc. Default is"euclidean".- calc_ENTRlogical, optional
A logical value indicating whether to calculate entropy (ENTR) along with DET. Default is
TRUE.- platformcharacter string
A character string specifying the computing platform. Options may include
"cpu","nv_gpu", etc.
- Returns:
- data.frame
A data frame containing:
Delay: Specific time delay from the values set in the parameters.
Embedding: Specific embedding dimension from the values set in the parameters.
Lmin: Minimal diagonal line lengths set for DET computation.
DET: The determinism values computed for the given input parameters.
ENTR (if
calc_ENTR = TRUE): The entropy values corresponding to the DET computations.RR: RR values.
- Example:
# Example usage
input_data <- matrix(runif(100), nrow = 10)
tau <- c(1, 2)
emb <- c(2, 3)
lmin <- 1
threshold <- 1
norm_method <- "euclidean"
calculate_entropy <- TRUE
comp_platform <- "cpu"
results <- accrqa_DET(
input = input_data,
tau_values = tau,
emb_values = emb,
lmin_values = lmin,
threshold_values = threshold,
norm = norm_method,
calc_ENTR = calculate_entropy,
platform = comp_platform
)
- accrqa_LAM(input, tau_values, emb_values, vmin_values, threshold_values, norm = "euclidean", calc_ENTR = TRUE, platform)
This function computes laminarity (LAM) based on the given input time series and RQA parameters.
Laminarity (LAM) is a measure in recurrence quantification analysis that describes the tendency of points to form vertical lines in the recurrence plot. This function provides configurable parameters for calculating LAM with options for normalization and entropy computation.
- Parameters:
- inputnumeric vector
A numeric vector representing the input time series data.
- tau_valuesnumeric vector
A numeric vector of time delay values.
- emb_valuesnumeric vector
A numeric vector of embedding dimension values.
- vmin_valuesnumeric vector
A numeric vector of minimum vertical line lengths.
- threshold_valuesnumeric vector
A numeric vector of threshold values for recurrence detection.
- normcharacter string, optional
A character string specifying the distance norm to use. Possible values are:
"euclidean": Euclidean distance."maximal": Maximum norm (Chebyshev distance)."none": No normalization.
Default is
"euclidean".- calc_ENTRlogical, optional
A logical value indicating whether to calculate entropy (
TRUEorFALSE). Default isTRUE.- platformcharacter string
A character string specifying the computing platform. Options may include
"cpu","nv_gpu", etc.
- Returns:
- data.frame
A data frame with the following columns:
LAM: Laminarity percentage.
V: Mean vertical line length.
Vmax: Maximum vertical line length.
ENTR: Entropy of the vertical line length distribution (if
calc_ENTR = TRUE).
- Example:
# Example usage of accrqa_LAM
input <- c(1.0, 2.0, 3.0, 4.0)
tau_values <- c(1, 2)
emb_values <- c(2, 3)
vmin_values <- c(2, 3)
threshold_values <- c(0.1, 0.2)
norm <- "euclidean"
calc_ENTR <- TRUE
result <- accrqa_LAM(
input,
tau_values,
emb_values,
vmin_values,
threshold_values,
norm,
calc_ENTR
)
- accrqa_RR(input, tau_values, emb_values, threshold_values, norm = "euclidean", platform)
This function computes the recurrence rate (RR) for a given input time series based on the specified delays, embedding dimensions, and thresholds. The function allows the user to specify normalization and computational platform.
Recurrence rate (RR) quantifies the density of recurrence points in a recurrence plot. This function uses a compiled C backend to efficiently compute RR based on the input parameters. It performs validations on input lengths and ensures that parameters like delays and embeddings are integers.
- Parameters:
- inputnumeric vector
A numeric vector representing the input time series.
- tau_valuesnumeric vector
A numeric vector of time delay (τ) values.
- emb_valuesnumeric vector
A numeric vector of embedding dimension (m) values.
- threshold_valuesnumeric vector
A numeric vector of threshold values for recurrence detection.
- normcharacter string, optional
A character string specifying the normalization method. Defaults to
"euclidean". Possible values are:"euclidean": Euclidean distance."maximal": Maximum norm (Chebyshev distance).
- platformcharacter string
A character string specifying the computational platform. Possible values are:
"cpu": Use the CPU for computations."nv_gpu": Use an NVIDIA GPU for computations.
- Returns:
- data.frame
A data frame containing:
Delay: The delay (τ) values used in the computation.
Embedding: The embedding dimension (m) values used.
Threshold: The threshold values used.
RR: The computed recurrence rate for each combination of parameters.
- Example:
# Example usage of accrqa_RR
input <- c(1.0, 2.0, 3.0, 4.0)
tau_values <- c(1, 2)
emb_values <- c(2, 3)
threshold_values <- c(0.1, 0.2)
norm <- "euclidean"
platform <- "cpu"
result <- accrqa_RR(
input,
tau_values,
emb_values,
threshold_values,
norm,
platform
)
print(result)