Model metrics#
- class sails.mvar_metrics.FourierMvarMetrics[source]#
This class for computing various connectivity metrics based on the Fourier transform of the MVAR co-efficients.
This is typically called using the
initialise()
classmethod to compute the frequency transform and return a class instance with methods for computing different connectivity estimators- classmethod initialise(model, sample_rate, freq_vect, nmodes=None)[source]#
Compute some basic values from the Fourier transform of the A matrix. This is a class method which creates and returns a
FourierMvarMetrics
instance with methods computing various connectivity metrics.This class computes the transfer function H the Fourier transform.
- Parameters:
model (sails LinearModel) – SAILS class instance containing a fitted linear model.
sample_rate (float) – Sample rate of the dataset
freq_vect (ndarray) – Vector of specifying which frequencies to compute connectivity at.
- Return type:
sails.FourierMvarMetrics instance
- class sails.mvar_metrics.ModalMvarMetrics[source]#
This class for computing various connectivity metrics based on a Modal decomposition of a fitted MVAR model
This is typically called using the
initialise()
orinitialise_from_modes()
classmethods to compute the decomposition and return a class instance containing them.- classmethod initialise(model, sample_rate, freq_vect, sum_modes=True)[source]#
Compute some basic values from the pole-residue decomposition of the A matrix. This is a class method which creates and returns a
ModalMvarMetrics
instance containing the decomposition parameters with methods computing various connectivity metrics.This class computes the transfer function H using modal parameters.
Currently only implemented for a single realisation (A.ndim == 3)
- Parameters:
model (sails LinearModel) – SAILS class instance containing a fitted linear model.
sample_rate (float) – Sample rate of the dataset
freq_vect (ndarray) – Vector of specifying which frequencies to compute connectivity at.
sum_modes (bool) – Flag indicating whether to sum across modes (Default value = True)
- Return type:
sails.ModalMvarMetrics instance
- classmethod initialise_from_modes(model, modes, sample_rate, freq_vect, mode_inds=None, sum_modes=True)[source]#
Compute some basic values from the pole-residue decomposition of the A matrix. This is a class method which creates and returns a
sails.ModalMvarMetrics
instance containing the decomposition parameters with methods computing various connectivity metrics.This class computes the transfer function H using modal parameters. A sub-set of modes can be used to compute a reduced H.
Currently only implemented for a single realisation (A.ndim == 3)
- Parameters:
model (sails LinearModel) – SAILS class instance containing a fitted linear model.
modes (sails MvarModalDecomposition) – SAILS class instance containing a fitted modal decomposition
sample_rate (float) – The sample rate of the fitted data
freq_vect (ndarray) – Vector of specifying which frequencies to compute connectivity at.
mode_inds (ndarray) – Which modes to use when constructing H (Default value = None)
sum_modes (Boolean) – Flag indicating whether to sum across modes when constructing H (Default value = True)
- Return type:
sails.ModalMvarMetrics instance
- sails.mvar_metrics.modal_transfer_function(evals, evecl, evecr, nchannels, sample_rate=None, freq_vect=None)[source]#
Compute the transfer function in pole-residue form, splitting the system into modes with separate transfer functions. The full system transfer function is then a linear sum of each modal transfer function.
- Parameters:
evals (ndarray) – Complex valued eigenvalues from an eigenvalue decomposition of an MVAR parameter matrix.
evecl (ndarray) – Complex valued left-eigenvectorsfrom an eigenvalue decomposition of an MVAR parameter matrix.
evecr (ndarry) – Complex valued right-eigenvectorsfrom an eigenvalue decomposition of an MVAR parameter matrix.
nchannels (int) – Number of channels in the decomposed system
sample_rate (float) – The sampling rate of the decomposed system (Default value = None)
freq_vect (ndarray) – Vector of frequencies at which to evaluate the transfer function (Default value = None)
- Returns:
Modal transfer-function (H) of size [nchannels x nchannels x nfrequencies x nmodes]
- Return type:
ndarray
Metric Mathemetical Functions#
- sails.mvar_metrics.sdf_spectrum(A, sigma, sample_rate, freq_vect)[source]#
Estimate of the Spectral Density as found on wikipedia and [Quirk1983].
This assumes that the spectral representation of A is invertable
- Parameters:
A (ndarray) – Matrix of autoregressive parameters, of size [nchannels x nchannels x model order]
sigma (ndarray) – Residual covariance matrix of the modelled system
sample_rate (float) – The samplingfrequency of the modelled system
freq_vect (ndarray) – Vector of frequencies at which to evaluate the spectrum
- Returns:
Frequenecy transform of input A, of size [nchannels x nchannels x nfrequencies]
- Return type:
ndarray
- sails.mvar_metrics.psd_spectrum(A, sigma, sample_rate, freq_vect)[source]#
Estimate the PSD representation of a set of MVAR coefficients as stated in [Penny2009] section 7.4.4.
This assumes that the spectral representation of A is invertable
WARNING: does not behave as expected for some data, use with caution. ar_spectrum is recommended.
- Parameters:
A (ndarray) – Matrix of autoregressive parameters, of size [nchannels x nchannels x model order]
sigma (ndarray) – Residual covariance matrix of the modelled system
sample_rate (float) – The samplingfrequency of the modelled system
freq_vect (ndarray) – Vector of frequencies at which to evaluate the spectrum
- Returns:
Frequenecy transform of input A, of size [nchannels x nchannels x nfrequencies]
- Return type:
ndarray
- sails.mvar_metrics.ar_spectrum(A, sigma, sample_rate, freq_vect)[source]#
Estimate the spectral representation of a set of MVAR coefficients as suggested by [Baccala2001], the equation without a number just below equation 13.
- Parameters:
A (ndarray) – Matrix of autoregressive parameters, of size [nchannels x nchannels x model order]
sigma (ndarray) – Residual covariance matrix of the modelled system
sample_rate (float) – The samplingfrequency of the modelled system
freq_vect (ndarray) – Vector of frequencies at which to evaluate the spectrum
- Returns:
Frequenecy transform of input A, of size [nchannels x nchannels x nfrequencies]
- Return type:
ndarray
- sails.mvar_metrics.transfer_function(Af)[source]#
Function for computing the transfer function of a system from the frequency transform of the autoregressive parameters.
- Parameters:
Af – Frequency domain version of parameters (can be calculated using ar_spectrum function) [nchannels x nchannels x nfreqs x nepochs]
- Returns:
System transfer function, of size [nchannels x nchannels x nfreqs x nepochs]
- Return type:
ndarray
Note
This function computes this the transfer function \(H\) as defined by:
\[H(f) = ( I - Af(f) )^{-1}\]where \(Af\) is the frequency transformed autoregessive parameter matrix and \(I\) is identity.
- sails.mvar_metrics.spectral_matrix(H, noise_cov)[source]#
Function for computing the spectral matrix, the matrix of spectra and cross-spectra.
- Parameters:
H (ndarray) – The transfer matrix [nchannels x nchannels x nfreqs x nepochs]
noise_cov (ndarray) – The noise covariance matrix of the system [nchannels x nchannels x nepochs]
- Returns:
System spectral matrix, of size [nchannels x nchannels x nfreqs x nepochs]
- Return type:
ndarray
Note
This function computes the system spectral matrix as defined by:
\[S(f) = H(f)\Sigma H(f)^H\]where \(H\) is the transfer function, \(\Sigma\) is the residual covariance and \(^H\) the Hermitian transpose.
- sails.mvar_metrics.coherency(S)[source]#
Method for computing the Coherency. This is the complex form of coherence, from which metrics such as magnitude squared coherence can be derived
- Parameters:
S (ndarray) – The system spectral matrix in 3D or 4D form, of size [nchannels x nchannels x nfreqs x nepochs]
- Returns:
Complex-valued coherency matrix, of size [nchannels x nchannels x nfreqs x nepochs]
- Return type:
ndarray
Note
This function computes the coherency as defined by:
\[Coh_{ij}(f) = \frac{S_{ij}(f)}{ \sqrt{ | S_{ii}(f)S_{jj}(f) | } }\]where \(S\) is the system spectral matrix.
- sails.mvar_metrics.partial_coherence(inv_S)[source]#
Method for computing the Partial Coherence.
- Parameters:
inv_S (ndarray) – The system spectral matrix in 3D or 4D form, of size [nchannels x nchannels x nfreqs x nepochs]
- Returns:
Complex-valued coherency matrix, of size [nchannels x nchannels x nfreqs x nepochs]
- Return type:
ndarray
Note
This function computes the partial coherence using the inverse spectral matrix \(P = S^{-1}\) where \(S\) is the system spectral matrix. The partial coherence is then computed as:
\[PCoh_{ij}(f) = \frac{ |P_{ij}(f)|^2 }{ \sqrt{ | P_{ii}(f)P_{jj}(f) | } }\]
- sails.mvar_metrics.partial_directed_coherence(Af)[source]#
Function to estimate the partial directed coherence from a set of multivariate parameters as defined in [Baccala2001].
- Parameters:
Af (ndarray) – Frequency domain version of parameters (can be calculated using ar_spectrum function), of size [nchannels x nchannels x nfrequencies x nrealisations]
- Returns:
The partial directed coherence of the system.
- Return type:
ndarray
Note
This function computes the partial coherence using the frequency transform of the autoregressive parameters \(Af\) subtracted from Identity.
\[\bar{Af}(f) = I - Af(f)\]The partial directed coherence is then
\[PDC_{ij}(f) = \frac{ | \bar{Af}_{ij}(f)| }{ \sqrt{ \sum_i | \bar{Af}_{ij}(f) |^2 } }\]
- sails.mvar_metrics.isolated_effective_coherence(Af, noise_cov)[source]#
Function for estimating the Isolated Effective Coherence as defined in [Pascual-Marqui2014].
- Parameters:
Af (ndarray) – Frequency domain version of parameters (can be calculated using ar_spectrum function) [nchannels x nchannels x nfreqs x nepochs]
noise_cov (ndarray) – The noise covariance matrix of the system [nchannels x nchannels x nepochs]
- Returns:
The isolated effective coherence.
- Return type:
ndarray
- sails.mvar_metrics.directed_transfer_function(H)[source]#
Method for computing the Directed Transfer Function as defined in [Kaminski1991].
- Parameters:
H (ndarray) – The transfer matrix [nchannels x nchannels x nfreqs x nepochs]
- Returns:
The directed transfer function of the system.
- Return type:
ndarray
Note
This function computes the directed transfer function from the system transfer function \(H\).
\[DTF_{ij}(f) = \frac{ | H_{ij}(f) |^2 }{ \sqrt{ \sum_j | \bar{H}_{ij}(f) |^2 } }\]
- sails.mvar_metrics.geweke_granger_causality(S, H, sigma)[source]#
This function computes the Geweke-Granger causality as defined in [Barrett2010]
- Parameters:
S (ndarray) – Spectral matrix [nchannels x nchannels x nfreqs x nepochs]
H (ndarray) – The transfer matrix [nchannels x nchannels x nfreqs x nepochs]
sigma (ndarray) – Residual noise covariance matrix
- Returns:
The Geweke-Granger causality
- Return type:
ndarray