cedalion.sigdecomp.multimodal.tcca_models
Module for temporally embedded CCA-like models. The temporally embedded technique is based on Bießmann et al. [BMG+10]
Functions
|
Construct a time-embedded version of a matrix X. |
Classes
|
Perform Elastic Net Canonical Correlation Analysis (CCA) between X_emb and Y. |
Class for decomposing multimodal data (X and Y) into latent sources using linear filters, using temporal embedding. |
|
|
Perform structured sparse Canonical Correlation Analysis (ssCCA) between two datasets X_emb and Y. |
|
Perform tCCA between two datasets X and Y. |
- class cedalion.sigdecomp.multimodal.tcca_models.MultimodalSourceDecompositionWithTemporalEmbedding(
- N_components: int = None,
- time_shifts: ndarray = None,
- max_iter: int = 100,
- tol: float = 1e-06,
- scale: bool = True,
- shift_source=True,
Bases:
MultimodalSourceDecomposition
Class for decomposing multimodal data (X and Y) into latent sources using linear filters, using temporal embedding.
This main class is inherited by other source decomposition methods, such as tCCA. It implements methods to validate input dimensions, apply normalization, and transform data from two modalities using filters learned during training. It assumes modality Y is delayed with respect to modality X.
- Parameters:
N_components (int, optional) – Number of components to extract. If None, the number of components is set to the minimum number of features between modalities.
time_shifts (np.ndarray) – Array with time shifts to be used for temporal embedding.
max_iter (int) – Maximum number of iterations for the algorithm.
tol (float) – Tolerance for convergence.
scale (bool) – Whether to scale the data during normalization to unit variance. Defaults to True.
shift_source (bool) – Whether to shift the reconstructed sources by the optimal time shift found during training. Defaults to True.
- estimate_optimal_shift(
- X_emb: DataArray,
- Y: DataArray,
Find optimal time shifts for X, one per component.
It finds the optimal time shift for each component by looking for the time-shifted X that produces the biggest correlation between reconstructed sources sx and sy.
- Parameters:
X_emb (DataArray) – Time-embedded version of X with dimensions (time_shift, sample_name, featureX_name).
Y (DataArray) – Input data for modality Y with dimensions (sample_name, featureY_name).
- shift_by_optimal(
- X: DataArray,
Shift X by optimal time shift using zero padding.
- transform(
- X: DataArray,
- Y: DataArray,
Apply the linear transformation on the input data using learnt filters.
This method validates the dimension labels and sizes of the input data to ensure consistency with the training data, perform temporal embedding on X, applies normalization using the stored parametersand, and then projects the normalized data onto a lower-dimensional space using the learned filters Wx and Wy. It retrieves the transformed arrays, a.k.a reconstructed sources.
- Parameters:
X (DataArray) – Input data for modality X. Expected to have dimensions (sample_name, featureX_name).
Y (DataArray) – Input data for modality Y. Expected to have dimensions (sample_name, featureY_name).
- Returns:
- A tuple (X_new, Y_new) where:
X_new (DataArray): Transformed data for modality X. Y_new (DataArray): Transformed data for modality Y.
- Return type:
tuple
- class cedalion.sigdecomp.multimodal.tcca_models.ElasticNetTCCA(
- N_components: int = None,
- l1_reg: float | list[float, float] = 0,
- l2_reg: float | list[float, float] = 0,
- time_shifts=None,
- max_iter: int = 100,
- tol: float = 1e-06,
- scale: bool = True,
- pls: bool = False,
- shift_source=True,
Bases:
MultimodalSourceDecompositionWithTemporalEmbedding
Perform Elastic Net Canonical Correlation Analysis (CCA) between X_emb and Y.
Apply temporally embedded CCA (tCCA) with L1 + L2 regularization, a.k.a elastic net. The algorithm finds sparse (L1) and normalized (L2) vectors Wx, and Wy as the solution to the following constrained optimization problem:
maximize Wx^T Cxy Wy subject to Wx^T Cx Wx = 1, Wy^T Cy Wy = 1,
||Wx||_1 <= c1x, ||Wy||_1 <= c1y, ||Wx||^2_2 <= c2x, ||Wy||^2_2 <= c2y
where Cx, Cy, and Cxy are the individual and cross-covariance matrices between X_emb and Y datasets, and the last four constraints correspond to the standard L1-norm and L2-norm penalization terms. c1x and c1y controls sparsity while c2x and c2y controls the magnitude of the vectors. PLS algorithms are also captured by this algorithm by sending Cx and Cy to the identity matrices.
The temporally embedded matrix X_emb is constructed by concatenating time-shifted versions of the original X. Y is assumed to be the delayed signals with respect to X so time shifts are always positive.
For the one-unit algorithm (sparse) SVD decomposition is performed on the whitened cross-covariance matrix K = Cx^(-1/2) Cxy Cy^(-1/2) (reduced to K = Cxy for PLS), using the following standard alternating power method (based on Parkhomenko et al. [PTB09]):
- Update u:
u <- K * v
u <- u / ||u||
- If L1:
u <- SoftThresholding(u, lambda_u/2) u <- u / ||u||
- Update v:
v <- K^T * u
v <- v / ||v||
- If L1:
v <- SoftThresholding(v, lambda_v/2) v <- v / ||v||
The resulting u and v are the leading left and right singular vectors of K which are nothing but individual components of the filters Wx and Wy. The softthresholding function bring some components to zero. If L2 regularization is used, prior to computing K, Cx and Cy are shifted by Cx <- Cx + alpha_x I and Cy <- Cy + alpha_y I.
Multiple components are obtained via a deflation method, subtracting from K its 1-rank approximation on each iteration. The returned vectors Wx and Wy are ordered in desceding order w.r.t. the singular values, which coincide with the canonical correlations.
- Parameters:
N_components (int, optional) – Number of components to extract. If None, the number of components is set to the minimum number of features between modalities.
l1_reg (float or list of floats) – list containing lambda_u and lambda_v (see above). If a single float is provided,
lambda_v. (then lambda_u =)
l2_reg (float or list of floats) – list containing alpha_x and alpha_y (see above). If a single float is provided,
alpha_y. (then alpha_x =)
time_shifts (np.ndarray) – Array with time shifts to be used for temporal embedding.
max_iter (int) – Maximum number of iterations for the algorithm.
tol (float) – Tolerance for convergence.
scale (bool) – Whether to scale the data during normalization to unit variance. Defaults to True.
pls (bool) – Whether to perform PLS regression. Defaults to False.
shift_source (bool) – Whether to shift the reconstructed sources by the optimal time shift found during training.
- Wx[source]
Linear filters for dataset X with dimensions (featureX_name, latent_featureX_name)
- Type:
DataArray
- Wy[source]
Linear filters for dataset Y with dimensions (featureY_name, latent_featureY_name).
- Type:
DataArray
- fit(
- X: DataArray,
- Y: DataArray,
- featureX_name: str = 'channel',
- featureY_name: str = 'channel',
Find the canonical vectors Wx, and Wy for the datasets X and Y.
- Parameters:
X (DataArray) – Input data for modality X. Expected to have dimensions (sample_name, featureX_name).
Y (DataArray) – Input data for modality Y. Expected to have dimensions (sample_name, featureY_name).
featureX_name (str, optional) – Label for X-feature dimension, set to ‘channel’ by default.
featureY_name (str, optional) – Label for Y-feature dimension, set to ‘channel’ by default.
- class cedalion.sigdecomp.multimodal.tcca_models.StructuredSparseTCCA(
- N_components: int = None,
- Lx=None,
- Ly=None,
- time_shifts=None,
- l1_reg: float | list[float, float] = 0,
- l2_reg: float | list[float, float] = 0,
- max_iter: int = 100,
- tol: float = 1e-06,
- scale: bool = True,
- shift_source: bool = True,
- pls: bool = False,
Bases:
MultimodalSourceDecompositionWithTemporalEmbedding
Perform structured sparse Canonical Correlation Analysis (ssCCA) between two datasets X_emb and Y.
The sstCCA algorithm is a temporally embedded extension of Chen et al. [CBL+13], and it assumes the underlying X and Y features are linked through a graph structure. It finds sparse (L1) vectors Wx, and Wy as the solution to the following constrained optimization problem:
maximize Wx^T Cxy Wy subject to Wx^T Cx Wx = 1, Wy^T Cy Wy = 1,
||Wx||_1 <= c1x, ||Wy||_1 <= c1y, Wx^T Lx Wx <= c2x, Wy^T Ly Wy <= c2y
where Cx, Cy, and Cxy are the individual and cross-covariance matrices between X_emb and Y datasets. The second constraint is the standard L1-norm penalization term, while the last constraint incorporates local information of the spatial distribution of the features trough the Laplacian matrices Lx and Ly. These terms encaurage filter components that are linked on the graphical structure to have similar values, making them to vary smoothly across the graph. The c1x and c1y controls sparsity while c2x and c2y controls the relative importante of the graph structure.
The temporally embedded matrix X_emb is constructed by concatenating time-shifted versions of the original X. Y is assumed to be the delayed signals with respect to X so time shifts are always positive.
For the one-unit algorithm, first Cx and Cy are shifted by Cx <- Cx + alpha_x Lx and Cy <- Cy + alpha_y Ly, and then SVD decomposition is performed on the whitened cross-covariance matrix K = Cx^(-1/2) Cxy Cy^(-1/2), using the following standard alternating power method (based on Parkhomenko et al. [PTB09]):
- Update u:
u <- K * v
u <- u / ||u||
- If L1:
u <- SoftThresholding(u, lambda_u/2) u <- u / ||u||
- Update v:
v <- K^T * u
v <- v / ||v||
- If L1:
v <- SoftThresholding(v, lambda_v/2) v <- v / ||v||
The resulting u and v are the leading left and right singular vectors of K which are nothing but individual components of the filters Wx and Wy. The softthresholding function bring some components to zero.
Multiple components are obtained via a deflation method, subtracting from K its 1-rank approximation on each iteration. The returned vectors Wx and Wy are ordered in desceding order w.r.t. the singular values, which coincide with the canonical correlations.
- Parameters:
N_components (int, optional) – Number of components to extract. If None, the number of components is set to the minimum number of features between modalities.
Lx (ndarray) – Laplacian matrix for modality X.
Ly (ndarray) – Laplacian matrix for modality Y.
time_shifts (np.ndarray) – Array with time shifts to be used for temporal embedding.
l1_reg (float or list of floats) – list containing lambda_u and lambda_v (see above). If a single float is provided,
lambda_v. (then lambda_u =)
l2_reg (float or list of floats) – list containing alpha_x and alpha_y (see above). If a single float is provided,
alpha_y. (then alpha_x =)
max_iter (int) – Maximum number of iterations for the algorithm.
tol (float) – Tolerance for convergence.
scale (bool) – Whether to scale the data during normalization to unit variance. Defaults to True.
pls (bool) – Whether to perform PLS regression. Defaults to False.
shift_source (bool) – Whether to shift the reconstructed sources by the optimal time shift found during training.
- Wx[source]
Linear filters for dataset X with dimensions (featureX_name, latent_featureX_name)
- Type:
DataArray
- Wy[source]
Linear filters for dataset Y with dimensions (featureY_name, latent_featureY_name).
- Type:
DataArray
- fit(
- X: DataArray,
- Y: DataArray,
- featureX_name: str = 'channel',
- featureY_name: str = 'channel',
Find the canonical vectors Wx, and Wy for the datasets X and Y.
- Parameters:
X (DataArray) – Input data for modality X. Expected to have dimensions (sample_name, featureX_name).
Y (DataArray) – Input data for modality Y. Expected to have dimensions (sample_name, featureY_name).
featureX_name (str, optional) – Label for X-feature dimension, set to ‘channel’ by default.
featureY_name (str, optional) – Label for Y-feature dimension, set to ‘channel’ by default.
- class cedalion.sigdecomp.multimodal.tcca_models.tCCA(
- N_components: int = None,
- time_shifts: ndarray = None,
- max_iter: int = 100,
- tol: float = 1e-06,
- scale: bool = True,
- shift_source=True,
Bases:
ElasticNetTCCA
Perform tCCA between two datasets X and Y.
This algorithm is a particular case of the one implemented in the ElasticNetTCCA class. See there for a detailed explanation of the algorithm.
- Parameters:
N_components (int, optional) – Number of components to extract. If None, the number of components is set to the minimum number of features between modalities.
time_shifts (np.ndarray) – Array with time shifts to be used for temporal embedding.
max_iter (int) – Maximum number of iterations for the algorithm.
tol (float) – Tolerance for convergence.
scale (bool) – Whether to scale the data during normalization to unit variance. Defaults to True.
shift_source (bool) – Whether to shift the reconstructed sources by the optimal time shift found during training.
- Wx[source]
Linear filters for dataset X with dimensions (featureX_name, latent_featureX_name)
- Type:
DataArray
- cedalion.sigdecomp.multimodal.tcca_models.temporal_embedding(
- X: DataArray,
- time_shifts: ndarray,
Construct a time-embedded version of a matrix X.
If X has shape T x N, the embedding matrix X_emb has shape P x T x N, where P is the number of time shifts in time_shifts. The embedding matrix is built by concatenating time-shift versions of the original matrix using the shifts inside time_shifts. Zero padding is used at the beggining of each time-shifted copy to preserve the original length of time direction.
- Parameters:
X (DataArray) – Input data with at least a time dimension.
time_shifts (np.ndarray) – Array with time shifts to be used for temporal embedding.
- Returns:
Time-embedded version of X with dimensions (time_shift, time, …).
- Return type:
DataArray