Synchrony analysis tools

Functions:

FCD_from_envelopes(X[, f_low, f_high, fs, ...])

Extract the FCD from the low-pass band filtered envelopes

FC_filtered(X[, f_low, f_high, fs, ...])

Calculates the Functional Connectivity matrix from the low-pass filterd envelopes of the signals in X.

FC_filtered_windowed(X[, t_start, t_end, ...])

Calculates the Functional Connectivity matrix from the low-pass filterd envelopes of the signals in X.

KuramotoOrderParameter(X, δ=None, σ=None[, ...])

Kuramoto order parameter, optional: subnetwork with indexes δ optional: frequency correction σ=stim frequency, T: tmax, dt:Ts

absDiffPhase(x)

The absolute difference between phases time series.

coherence(x[, nperseg, noverlap, fs, applySin])

Coherence for each pair of signals in x

completeSynchrony(x1, x2)

Absolute error between two signals.

complex_coherence(data_x, data_y[, nfft, ...])

Squared complex coherence, from here is easy to obtain the absolute, the real or the imaginary value

complex_coherence_matrix(data[, nfft, ...])

Squared complex coherence, from here is easy to obtain the absolute, the real or the imaginary value

conditionalProbabiliy(x1, x2[, theta, ...])

Conditional probability

cross_spectrum(x[, nperseg, noverlap, fs, ...])

Real part of the cross-spectrum for each pair of signals in x

diffPhaseHilbert(X[, f_low, f_high, fs, ...])

Difference of phases from applying the Hilbert Transform in the signals of X.

diffPhaseTheta(X)

Difference of phases assuming that X contains phase time series.

durationfromLabels(labels[, time_window, ...])

Calculate the duration of the events from a list of the labels assigned to each time window.

elementFunctionalConnectivity(x1, x2)

Calculate the element i,j of the functional connectivity matrix

elementFunctionalConnectivityTheta(x1, x2)

Calculate the i,j element of the functional connectivity matrix after applying the sin() function to the data's vectors.

entropySynchrony(x1[, x2, n, m])

Index of synchrony based in the Shannon entropy If bin_start and bin_end, are not specified, this functions actuates similar as np.histogram to get the histogram of the i,j element of the matrix X.

entropySynchronyMatrix(X[, nbins, ...])

Calculate the entropy index of synchrony

envelopesFrequencyBand(X[, f_low, f_high, ...])

Returns the envelpes of the Hilbert transform of the signal at a specific frequency

extractEventsTimes(x[, high_value, min_duration])

Returns a list with the starting and ending time points from the binary time serie of the thresholded events

extractTimeStatisticsEvents(X[, min_duration])

Extract the co-occurrrences of events in several nodes And their time characteristics: fractional occupancy in each node and durations.

extract_FCD(data[, wwidth, maxNwindows, ...])

Created on Wed Apr 27 15:57:38 2016 @author: jmaidana @author: porio

fourierModeIndex(x1, x2[, n, m])

Fourier index based in the trigonometry identity sin^2(alpha)+cos^2 (alpha)=1

high_order_cooccurrences(events_times)

Find the co-occurrences

hilbertFrequencyBand(X[, f_low, f_high, fs, ...])

Calculate to envelope and phase of a signal in the frequency band specified by [f_loww , f_high ] Hz.

localOrderParameter(x)

Local order parameter Relation of each node phase with the global order parameter

phaseLockingDiffPhase(diffX)

Phase locking value given the phase difference matrix

phaseLockingValueMatrix(X)

Phase locking value for pairs in matrix of phases X

phaseLockingValueTwoNodes(x1, x2)

Phase locking between x1 and X2

shannonEntropy(p)

Calculates the Shannon entropy of the probability mass distribution p[j]

shannonEntropyTimeMatrix(X[, nbins, ...])

Calculates the shanon entropy of the probability mass distribution p[j].

synchronyMatrices(X[, start_time, end_time])

Return multiple synchrony measurements for a data matrix with size NxT

synchronyTwoNodes(x1, x2[, n, m])

Return multiple synchrony measurements for two vector of data with the same length

transitionsfromLabels(labels)

Calculate the transitions between events from a list of the labels assigned to each time window.

analysis.synchronization.FCD_from_envelopes(X, f_low=8, f_high=13, fs=1000, wwidth=1000, olap=0.9, mode='corr', applyTrim=True, applySin=True)

Extract the FCD from the low-pass band filtered envelopes

Parameters
  • X (float 2D array) – Data matrix of size NxT.

  • f_low (float, optional) – Low frequency limit of the pass-band filter. The default is 5 Hz.

  • f_high (float, optional) – High frequency limit of the pass-band filter. The default is 13 Hz.

  • fs (int, optional) – Sampling frequency. The default is 1000 samples/second.

  • wwidth (integer, optional) – Length of data windows in which the series will be divided to calculate each FC, in samples the default is 1000

  • olap (float between 0 and 1, optional) – Overlap between neighboring data windows, in fraction of window length

  • mode (str, 'corr' | 'psync' | 'plock' | 'tdcorr') – Measure to calculate the Functional Connectivity (FC) between nodes. ‘corr’ : Pearson correlation. Uses the corrcoef function of numpy. ‘psync’ : Pair-wise phase synchrony. ‘plock’ : Pair-wise phase locking. ‘tdcorr’ : Time-delayed correlation, looks for the maximum value in a cross-correlation of the data series

  • applyTrim (boolean, optional) – Defines if the data comes from simulation, then the impulse response time is removed. The default is True.

  • applySin (boolean, optional) – Defines if the sin function must be applied to the data before any processing.

Returns

  • FCDmatrix (numpy array) – Correlation matrix between all the windowed FCs.

  • CorrVectors (numpy array) – Collection of FCs, linearized. Only the lower triangle values (excluding the diagonal) are returned

  • shift (integer) – The distance between windows that was actually used (in samples)

analysis.synchronization.FC_filtered(X, f_low=0.5, f_high=100, fs=1000, applyTrim=True, applySin=True)

Calculates the Functional Connectivity matrix from the low-pass filterd envelopes of the signals in X. The envelopes correspond to the frequency band specified by [f_low , f_high ] Hz.

Parameters
  • X (float 2D array) – Data matrix of size NxT.

  • f_low (float, optional) – Low frequency limit of the pass-band filter. The default is 0.5 Hz.

  • f_high (float, optional) – High frequency limit of the pass-band filter. The default is 100 Hz.

  • fs (int, optional) – Sampling frequency. The default is 1000 samples/second.

  • type (str, optional) – Type of the pass-band filter. The default is ‘butterworth’, the other option is ‘chebysev’.

  • applyTrim (boolean, optional) – Defines if the data comes from simulation, then the impulse response time is removed. The default is True.

  • applySin (boolean, optional) – Defines if the sin function must be applied to the data before any processing.

Returns

  • FC (float 2D array) – Functional connectivity matrix, uses the pearson coefficient, then the values are in the range [-1,1].

  • mean_energy (float) – Average of the energy from the envelopes of the signal in X.

analysis.synchronization.FC_filtered_windowed(X, t_start=20000, t_end=40000, f_low=0.5, f_high=100, fs=1000, applyTrim=True, applySin=True)

Calculates the Functional Connectivity matrix from the low-pass filterd envelopes of the signals in X. The envelopes correspond to the frequency band specified by [f_low , f_high ] Hz, and between the time points specified by [t_start , t_end ) samples.

Parameters
  • X (float 2D array) – Data matrix of size NxT.

  • f_low (float, optional) – Low frequency limit of the pass-band filter. The default is 0.5 Hz.

  • f_high (float, optional) – High frequency limit of the pass-band filter. The default is 100 Hz.

  • fs (int, optional) – Sampling frequency. The default is 1000 samples/second.

  • type (str, optional) – Type of the pass-band filter. The default is ‘butterworth’, the other option is ‘chebysev’.

  • applyTrim (boolean, optional) – Defines if the data comes from simulation, then the impulse response time is removed. The default is True.

  • applySin (boolean, optional) – Defines if the sin function must be applied to the data before any processing.

  • t_start (int, optional) – Initial time point. The default is 20000 samples.

  • t_end (int, optional) – Final time point. The default is 40000 samples.

Returns

  • FC (float 2D array) – Functional connectivity matrix, uses the pearson coefficient, then the values are in the range [-1,1].

  • mean_energy (float) – Average of the energy from the envelopes of the signal in X.

analysis.synchronization.KuramotoOrderParameter(X, δ=-1, σ=0, T=-1, dt=-1)

Kuramoto order parameter, optional: subnetwork with indexes δ optional: frequency correction σ=stim frequency, T: tmax, dt:Ts

Parameters
  • X (2D array) – Nodes x time.

  • δ (int 1D array, optional) – nodes indexes. The default is -1.

  • σ (float, optional) – stim frequency. The default is 0.

  • T (float, optional) – t_max. The default is -1.

  • dt (float, optional) – sampling time. The default is -1.

Returns

KOP – Kuramoto order parameter. Between 0 and 1.

Return type

float

analysis.synchronization.absDiffPhase(x)

The absolute difference between phases time series. The max value in the absolute difference between phases is pi.

Parameters

x (float (ND array)) – Array that contains difference of phases information.

Returns

Difference of phases in degrees. Useful for the visualization of synchrony with one node as reference.

Return type

abs_diff same as x

analysis.synchronization.coherence(x, nperseg=4096, noverlap=2048, fs=1000, applySin=False)

Coherence for each pair of signals in x

Parameters
  • x (float 2D array.) – Matrix with the signal of size Nodes x Time.

  • nperseg (int, optional) – Number of samples of the time window. The default is 4096.

  • noverlap (int, optional) – Number of samples of the overlap window. The default is 3600. With None, the nperseg is used as value of nfft of a FFT with zero-padding

  • fs (int, optional) – Sampling frequency

  • applySin (boolean, optional) – Apply the sin function before performing the coherence calculation

Returns

  • freqs (float 1D array) – frequency bins

  • Cxx (float 3D array) – Tensor of coherences of each pair of nodes. Node x Node x Frequency

analysis.synchronization.completeSynchrony(x1, x2)

Absolute error between two signals. From the use of native Python operators, the arrays could be of any dimension while np.shape ( x1 )==*np.shape* ( x2 )

Parameters
  • x1 (float 1D (ND) vector) – Array of data.

  • x2 (float 1D (ND) vector) – Array of data.

Returns

error – Absolute error.

Return type

float

analysis.synchronization.complex_coherence(data_x, data_y, nfft=5000, freq_index=1, wcoh=1000)

Squared complex coherence, from here is easy to obtain the absolute, the real or the imaginary value

Parameters
  • data_x (float array) – data from channel x

  • data_y (float array) – data from channel y

  • nfft (int, optional) – number of points of the FFT. Defines the resolution of the spectrums. The default is 5000.

  • freq_index (int or array int, optional) – indices of the frequencies of interest, relative to nfft. The default is 1.

  • wcoh (TYPE, optional) – time window. The default is 1000.

Returns

coh – Complex average of the coherence in the frequency of interest.

Return type

complex

analysis.synchronization.complex_coherence_matrix(data, nfft=5000, freq_index=1, wcoh=1000)

Squared complex coherence, from here is easy to obtain the absolute, the real or the imaginary value

Parameters
  • data_x (2D float array) – data with shape NxT

  • nfft (int, optional) – number of points of the FFT. Defines the resolution of the spectrums. The default is 5000.

  • freq_index (int or array int, optional) – indices of the frequencies of interest, relative to nfft. The default is 1.

  • wcoh (TYPE, optional) – time window. The default is 1000.

Returns

coh – Complex average of the coherence in the frequency of interest.

Return type

complex

analysis.synchronization.conditionalProbabiliy(x1, x2, theta=6.280043714525997, epsilon=0.001, n=1, m=1)

Conditional probability

Parameters
  • x1 (1D array) – float Lenghth L.

  • x2 (float 1D array) – Length L.

  • theta (float, optional) – Target-phase at where the signals must coincide. The default is 1.999*np.pi.

  • epsilon (float, optional) – error around the target-phase

  • n (int, optional) – Harmonic of x1 signal. The default is 1.

  • m (int, optional) – Harmonic of x2 signal. The default is 1.

Returns

eta – Percentage of coincident points

Return type

float

analysis.synchronization.cross_spectrum(x, nperseg=4096, noverlap=2048, fs=1000, applySin=False)

Real part of the cross-spectrum for each pair of signals in x

Parameters
  • x (float 2D array.) – Nodes x Time.

  • nperseg (int, optional) – Number of samples of the time window. The default is 4096.

  • noverlap (int, optional) – Number of samples of the overlap window. The default is 3600. With None, the nperseg is used as value of nfft of a FFT with zero-padding

  • fs (int, optional) – Sampling frequency

  • applySin (boolean, optional) – Apply the sin function before performing the cross-spectrum calculation

Returns

  • freqs (float 1D array) – frequency bins

  • Cxx (float 3D array) – Tensor of coherences of each pair of nodes. Node x Node x Frequency

analysis.synchronization.diffPhaseHilbert(X, f_low=2, f_high=100, fs=1000, applyTrim=True, applySin=False)

Difference of phases from applying the Hilbert Transform in the signals of X.

Parameters
  • X (float 2D array) – Data matrix of size NxT.

  • f_low (float, optional) – Low frequency limit of the pass-band filter. The default is 0.5 Hz.

  • f_high (float, optional) – High frequency limit of the pass-band filter. The default is 100 Hz.

  • fs (int, optional) – Sampling frequency. The default is 1000 samples/second.

  • type (str, optional) – Type of the pass-band filter. The default is ‘butterworth’, the other option is ‘chebysev’.

  • applyTrim (boolean, optional) – Defines if the data comes from simulation, then the impulse response time is removed. The default is True.

  • applySin (boolean, optional) – Defines if the sin function must be applied to the data before any processing.

Returns

  • angles float 2D array – Phase time series.

  • diff (float 2D arrat) – Difference of the phases. Size TxNxN

analysis.synchronization.diffPhaseTheta(X)

Difference of phases assuming that X contains phase time series.

Parameters

X (float 2D array) – Matrix of phase time series. Size N x T.

Returns

  • angles float 2D array – The transpose of X and applied a module of 2 pi

  • diff (float 2D arrat) – Difference of the phases. Size TxNxN

analysis.synchronization.durationfromLabels(labels, time_window=118, overlap=0.5)

Calculate the duration of the events from a list of the labels assigned to each time window. The method considers that the time windows could have overlap between them.

Parameters
  • labels (int array) – Array with the labels from 0 to N-1 of N groups, or clusters.

  • time_window (int, optional) – The length of the time window in time units. The default is 118 ms.

  • overlap (float, optional) – Overlap between the time windows. The default is 0.5 for 50% overlapping.

Returns

A dictionary which keys are the labels, and the items are the duration of the events of each label.

Return type

duration_clusters

analysis.synchronization.elementFunctionalConnectivity(x1, x2)

Calculate the element i,j of the functional connectivity matrix

Parameters
  • x1 (float 1D array) – First vector of data with lenghth L.

  • x2 (float 1D array) – Second vector of data with lenghth L.

Returns

fc – Pearson coefficient between x1 and x2, from -1 to 1.

Return type

float

analysis.synchronization.elementFunctionalConnectivityTheta(x1, x2)

Calculate the i,j element of the functional connectivity matrix after applying the sin() function to the data’s vectors.

Parameters
  • x1 (1D array: float) – First vector of data with lenghth L.

  • x2 (1D array: float) – Second vector of data with lenghth L.

Returns

fc – Pearson coefficient between x1 and x2, from -1 to 1.

Return type

float

analysis.synchronization.entropySynchrony(x1, x2=None, n=1, m=1)

Index of synchrony based in the Shannon entropy If bin_start and bin_end, are not specified, this functions actuates similar as np.histogram to get the histogram of the i,j element of the matrix X. :param x1: First vector of data with lenghth L. :type x1: float 1D array :param x2: First vector of data with lenghth L. :type x2: float 1D array. :param n: Harmonic of x1 signal. The default is 1. :type n: int, optional :param m: Harmonic of x2 signal. The default is 1. :type m: int, optional

Returns

rho – Synchrony index from 0 t0 1.

Return type

float

analysis.synchronization.entropySynchronyMatrix(X, nbins=100, bin_start=None, bin_end=None)

Calculate the entropy index of synchrony

Parameters
  • X (float 2D array) – matrix of data with size TxN.

  • nbins (int, optional) – Number of bins to discretize the element i,j of the tensor X. The default is 100 bins.

  • bin_start (float, optional) – Start value of the bins range. The default is None. If None

  • bin_end (float, optional) – End value of the bins range. The default is None.

Returns

rho – Entropy synchrony measurement, values from 0 to 1.

Return type

float

analysis.synchronization.envelopesFrequencyBand(X, f_low=0.5, f_high=100, fs=1000, applyTrim=True, applySin=True, applyLow=True, f_lowpass=0.5)

Returns the envelpes of the Hilbert transform of the signal at a specific frequency

Warning! envelopes has eight seconds lesser than X

Parameters
  • X (float 2D array) – Data matrix of size N x T.

  • f_low (float, optional) – Low frequency limit of the pass-band filter. The default is 0.5 Hz.

  • f_high (float, optional) – High frequency limit of the pass-band filter. The default is 100 Hz.

  • fs (int, optional) – Sampling frequency. The default is 1000 samples/second.

  • applyTrim (boolean, optional) – Defines if the data comes from simulation, then the impulse response time is removed. The default is True.

  • applyLow (boolean, optional) – Defines if the envelopes are low-pass filtered or if they are directly returned from the Hilbert transform.

  • f_lowpass (float, optional) – Defines the stop band of the low-pass filter applied to the Hilbert envelopes

  • applySin (boolean, optional) – Defines if the sin function must be applied to the data before any processing.

Returns

envelopes – Low-pass filtered envelopes. Size N x (T-5*fs)

Return type

2D array

analysis.synchronization.extractEventsTimes(x, high_value=1, min_duration=3)

Returns a list with the starting and ending time points from the binary time serie of the thresholded events

Parameters
  • x (1D int) – Thresholded events or excerpts of the signal at a specific frequency band.

  • high_value (int, optional) – The integer that codifies the event. The default is 1 for binary signals.

  • min_duration (int, optional) – The minimum duration of the event. The default is 3 samples.

Returns

Array of the events with the index of the node (row 0), starting time (row 1) and ending time (row 2). The size is 3 x N_events.

Return type

2D int array

analysis.synchronization.extractTimeStatisticsEvents(X, min_duration=5)

Extract the co-occurrrences of events in several nodes And their time characteristics: fractional occupancy in each node and durations.

Parameters
  • X (2D int array) – binary time series of size N(nodes) x T(sampling points). 1: indicate the presence of an event. 0: ausence

  • min_duration (float, optional) – Minimum overlap between the node’s events to be considered a co-occurrence. The default is 5 sample points.

Returns

  • durations (float) – Duration of each event.

  • occupancy (float) – Fractional occupancy of the events in each node N.

  • co_occurrences (list) – Each co-occurrence is a tuple with the starting time point, the duration, and a list with the nodes’ indexes.

analysis.synchronization.extract_FCD(data, wwidth=1000, maxNwindows=100, olap=0.9, nfft=5000, freq_index=1, wcoh=100, coldata=False, mode='corr')

Created on Wed Apr 27 15:57:38 2016 @author: jmaidana @author: porio

Source: https://github.com/vandal-uv/anarpy/blob/master/src/anarpy/utils/FCDutil/fcd.py Functional Connectivity Dynamics from a collection of time series

Parameters
  • data (array-like) – 2-D array of data, with time series in rows (unless coldata is True)

  • wwidth (integer) – Length of data windows in which the series will be divided, in samples

  • maxNwindows (integer) – Maximum number of windows to be used. wwidth will be increased if necessary

  • olap (float between 0 and 1) – Overlap between neighboring data windows, in fraction of window length

  • coldata (Boolean) – if True, the time series are arranged in columns and rows represent time

  • nfft (int) – Number of points of the FFT

  • freq_index (int or array int) – frequencies of interest

  • wcoh – internal window for the coherence

  • mode ('corr' | 'psync' | 'plock' | 'tdcorr') – Measure to calculate the Functional Connectivity (FC) between nodes. ‘corr’ : Pearson correlation. Uses the corrcoef function of numpy. ‘psync’ : Pair-wise phase synchrony. ‘plock’ : Pair-wise phase locking. ‘tdcorr’ : Time-delayed correlation, looks for the maximum value in a cross-correlation of the data series

Returns

  • FCDmatrix (numpy array) – Correlation matrix between all the windowed FCs.

  • CorrVectors (numpy array) – Collection of FCs, linearized. Only the lower triangle values (excluding the diagonal) are returned

  • shift (integer) – The distance between windows that was actually used (in samples)

analysis.synchronization.fourierModeIndex(x1, x2, n=1, m=1)

Fourier index based in the trigonometry identity sin^2(alpha)+cos^2 (alpha)=1

Parameters
  • x1 (1D array: float) – First vector of data with lenghth L.

  • x2 (1D array: float) – Second vector of data with lenghth L.

  • n (int, optional) – Harmonic of x1 signal. The default is 1.

  • m (int, optional) – Harmonic of x2 signal. The default is 1.

Returns

gamma

Return type

float.

analysis.synchronization.high_order_cooccurrences(events_times)

Find the co-occurrences

Parameters

events_times (2D int array) – Array of the events with the index of the node (row 0), starting time (row 1) and ending time (row 2). The size is 3 x N_events.

Returns

co_occurrences – List of the co-occurrent events. Each co-occurrence is a tuple with the starting time point, the duration, and a list with the nodes’ indexes.

Return type

list

analysis.synchronization.hilbertFrequencyBand(X, f_low=0.5, f_high=100, fs=1000, type='butterworth', applyTrim=True, applySin=True)

Calculate to envelope and phase of a signal in the frequency band specified by [f_loww , f_high ] Hz. Uses the hilbert transform, then the result has more accuraccy as narrower is the frequency band.

Parameters
  • X (float 2D array) – Data matrix of size NxT.

  • f_low (float, optional) – Low frequency limit of the pass-band filter. The default is 0.5 Hz.

  • f_high (float, optional) – High frequency limit of the pass-band filter. The default is 100 Hz.

  • fs (int, optional) – Sampling frequency. The default is 1000 samples/second.

  • type (str, optional) – Type of the pass-band filter. The default is ‘butterworth’, the other option is ‘chebysev’.

  • applyTrim (boolean, optional) – Defines if the data comes from simulation, then the impulse response time is removed. The default is True.

  • applySin (boolean, optional) – Defines if the sin function must be applied to the data before any processing.

Returns

  • amplitudes (1D array) – Envelope from the Hilbert transform. The length is 4*fs less than the length of X.

  • angles (1D array) – Phase time serie from the Hilbert transform. The length is 4*fs less than the length of X.

analysis.synchronization.localOrderParameter(x)

Local order parameter Relation of each node phase with the global order parameter

Parameters

x (2D array) – Nodes x time.

Returns

LOP – Array with values of the local order parameter.

Return type

float 1D array

analysis.synchronization.phaseLockingDiffPhase(diffX)

Phase locking value given the phase difference matrix

Parameters

diffX (float 3D array) – Time X Nodes X Nodes difference between phases

Returns

plv – phase locking value matrix

Return type

float 2D array

analysis.synchronization.phaseLockingValueMatrix(X)

Phase locking value for pairs in matrix of phases X

Parameters

X (float 2D array) – Array of data of size Nodes x Time

Returns

plv – Phase locking value matrix

Return type

float 2D array

analysis.synchronization.phaseLockingValueTwoNodes(x1, x2)

Phase locking between x1 and X2

Parameters
  • x1 (float 1D array) – First vector of data with lenghth L.

  • x2 (float 1D array) – Second vector of data with lenghth L.

Returns

plv – phase locking value

Return type

float

analysis.synchronization.shannonEntropy(p)

Calculates the Shannon entropy of the probability mass distribution p[j]

Parameters

p (float 1D array) – Probability mass distirbution. All the elements must sum 1.

Returns

S – Shannon Entropy.

Return type

float

analysis.synchronization.shannonEntropyTimeMatrix(X, nbins=100, bin_start=None, bin_end=None)

Calculates the shanon entropy of the probability mass distribution p[j]. If bin_start and bin_end, are not specified, this functions actuates similar as np.histogram to get the histogram of the i,j element of the tensor X.

Parameters
  • X (float 3D array) – A tensor of size TxMxN where T is the number of time points

  • nbins (int, optional) – Number of bins to discretize the element i,j of the tensor X. The default is 100 bins.

  • bin_start (float, optional) – Start value of the bins range. The default is None. If None

  • bin_end (float, optional) – End value of the bins range. The default is None.

Returns

S – Shannon Entropy matrix.

Return type

float 2D array

analysis.synchronization.synchronyMatrices(X, start_time=0, end_time=20000)

Return multiple synchrony measurements for a data matrix with size NxT

Parameters
  • X (float 2D array) – Data matrix of size N x T

  • start_time (int, optional) – Initial time point. The default is 0.

  • end_time (int, optional) – Final time point. The default is 20000.

Returns

  • plv_matrix (float 2D array) – Phase locking value array.

  • gamma_matrix (float 2D array) – Fourier mode index array.

  • phi_matrix (float 2D array) – Absolute error array.

  • SE_matrix (float 2D array) – Entropy Synchrony array.

analysis.synchronization.synchronyTwoNodes(x1, x2, n=1, m=1)

Return multiple synchrony measurements for two vector of data with the same length

Parameters
  • x1 (float 1D array) – First data vector with length L.

  • x2 (float 1D array) – Second data vector with length L.

  • n (int, optional) – Mode (Harmonic) of x1. The default is 1.

  • m (int, optional) – Mode (Harmonic) of x2. The default is 1.

Returns

  • plv (float) – Phase locking value.

  • gamma (float) – Fourier mode index.

  • phi (float) – Absolute error.

  • rho (float) – Entropy Synchrony.

analysis.synchronization.transitionsfromLabels(labels)

Calculate the transitions between events from a list of the labels assigned to each time window.

Parameters

labels (int array) – Array with the labels from 0 to N-1 of N groups, or clusters.

Returns

transition_matrix – Matrix with the counts of the transtions. The rows indicate the destine The columns indicate the orgin

Return type

2d int array