heat.decomposition
Add the decomposition functions to the ht.decomposition namespace
Submodules
Package Contents
- _isvd(new_data: heat.core.dndarray.DNDarray, U_old: heat.core.dndarray.DNDarray, S_old: heat.core.dndarray.DNDarray, Vt_old: heat.core.dndarray.DNDarray | None = None, maxrank: int | None = None, old_matrix_size: int | None = None, old_rowwise_mean: heat.core.dndarray.DNDarray | None = None) Tuple[heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray] | Tuple[heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray]
Helper function for iSVD and iPCA; follows roughly the “incremental PCA with mean update”, Fig.1 in: David A. Ross, Jongwoo Lim, Ruei-Sung Lin, Ming-Hsuan Yang. Incremental Learning for Robust Visual Tracking. IJCV, 2008.
Either incremental SVD / PCA or incremental SVD / PCA with mean subtraction is performed.
- Parameters:
new_data (DNDarray) – new data as DNDarray
U_old (DNDarrays) – “old” SVD-factors if no Vt_old is provided, only U and S are computed (PCA)
S_old (DNDarrays) – “old” SVD-factors if no Vt_old is provided, only U and S are computed (PCA)
Vt_old (DNDarrays) – “old” SVD-factors if no Vt_old is provided, only U and S are computed (PCA)
maxrank (int, optional) – rank to which new SVD should be truncated
old_matrix_size (int, optional) – size of the old matrix; this does not need to be identical to Vt_old.shape[1] as “old” SVD might have been truncated
old_rowwise_mean (int, optional) – row-wise mean of the old matrix; if not provided, no mean subtraction is performed
- class PCA(n_components: int | float | None = None, copy: bool = True, whiten: bool = False, svd_solver: str = 'hierarchical', tol: float | None = None, iterated_power: str | int = 0, n_oversamples: int = 10, power_iteration_normalizer: str = 'qr', random_state: int | None = None)[source]
Bases:
heat.TransformMixin,heat.BaseEstimatorPricipal Component Analysis (PCA).
Linear dimensionality reduction using Singular Value Decomposition of the data to project it to a lower dimensional space. The input data is centered but not scaled for each feature before applying the SVD.
- Parameters:
n_components (int, float, None, default=None) – Number of components to keep. If n_components is not set all components are kept. If n_components is an integer, it specifies the number of components to keep. If n_components is a float between 0 and 1, it specifies the fraction of variance explained by the components to keep.
copy (bool, default=True) – In-place operations are not yet supported. Please set copy=True.
whiten (bool, default=False) – Not yet supported.
svd_solver ({'full', 'hierarchical'}, default='hierarchical') – ‘full’ : Full SVD is performed. In general, this is more accurate, but also slower. So far, this is only supported for tall-skinny or short-fat data. ‘hierarchical’ : Hierarchical SVD, i.e., an algorithm for computing an approximate, truncated SVD, is performed. Only available for data split along axis no. 0. ‘randomized’ : Randomized SVD is performed.
tol (float, default=None) – Not yet necessary as iterative methods for PCA are not yet implemented.
iterated_power (int, default=0) – if svd_solver=’randomized’, this parameter is the number of iterations for the power method. Choosing iterated_power > 0 can lead to better results in the case of slowly decaying singular values but is computationally more expensive.
n_oversamples (int, default=10) – if svd_solver=’randomized’, this parameter is the number of additional random vectors to sample the range of X so that the range of X can be approximated more accurately.
power_iteration_normalizer ({'qr'}, default='qr') – if svd_solver=’randomized’, this parameter is the normalization form of the iterated power method. So far, only QR is supported.
random_state (int, default=None) – if svd_solver=’randomized’, this parameter allows to set the seed for the random number generator.
- Variables:
components (DNDarray of shape (n_components, n_features)) – Principal axes in feature space, representing the directions of maximum variance in the data. The components are sorted by explained_variance_.
explained_variance (DNDarray of shape (n_components,)) – The amount of variance explained by each of the selected components. Not supported by svd_solver=’hierarchical’ and svd_solver=’randomized’.
explained_variance_ratio (DNDarray of shape (n_components,)) – Percentage of variance explained by each of the selected components. Not supported by svd_solver=’hierarchical’ and svd_solver=’randomized’.
total_explained_variance_ratio (float) – The percentage of total variance explained by the selected components together. For svd_solver=’hierarchical’, an lower estimate for this quantity is provided; see
ht.linalg.hsvd_rtol()andht.linalg.hsvd_rank()for details. Not supported by svd_solver=’randomized’.singular_values (DNDarray of shape (n_components,)) – The singular values corresponding to each of the selected components. Not supported by svd_solver=’hierarchical’ and svd_solver=’randomized’.
mean (DNDarray of shape (n_features,)) – Per-feature empirical mean, estimated from the training set.
n_components (int) – The estimated number of components.
n_samples (int) – Number of samples in the training data.
noise_variance (float) – not yet implemented
Notes
Hierarchical SVD (svd_solver = “hierarchical”) computes an approximate, truncated SVD. Thus, the results are not exact, in general, unless n_components chosen is larger than the actual rank (=matrix rank) of the underlying data; see
ht.linalg.hsvd_rank()andht.linalg.hsvd_rtol()for details. Randomized SVD (svd_solver = “randomized”) is a stochastic algorithm that computes an approximate, truncated SVD.- n_components = None
- copy = True
- whiten = False
- svd_solver = 'hierarchical'
- tol = None
- iterated_power = 0
- n_oversamples = 10
- power_iteration_normalizer = 'qr'
- random_state = None
- components_ = None
- explained_variance_ = None
- explained_variance_ratio_ = None
- total_explained_variance_ratio_ = None
- singular_values_ = None
- mean_ = None
- n_components_ = None
- n_samples_ = None
- noise_variance_ = None
- fit(X: heat.DNDarray, y=None) Self[source]
Fit the PCA model with data X.
- Parameters:
X (DNDarray of shape (n_samples, n_features)) – Data set of which PCA has to be computed.
y (Ignored) – Not used, present for API consistency by convention.
- transform(X: heat.DNDarray) heat.DNDarray[source]
Apply dimensionality based on PCA to X.
- Parameters:
X (DNDarray of shape (n_samples, n_features)) – Data set to be transformed.
- inverse_transform(X: heat.DNDarray) heat.DNDarray[source]
Transform data back to its original space.
- Parameters:
X (DNDarray of shape (n_samples, n_components)) – Data set to be transformed back.
- class IncrementalPCA(n_components: int | None = None, copy: bool = True, whiten: bool = False, batch_size: int | None = None)[source]
Bases:
heat.TransformMixin,heat.BaseEstimatorIncremental Principal Component Analysis (PCA).
This class allows for incremental updates of the PCA model. This is especially useful for large data sets that do not fit into memory.
An example how to apply this class is given in, e.g., benchmarks/cb/decomposition.py.
- Parameters:
n_components (int, optional) – Number of components to keep. If n_components is not set all components are kept (default).
copy (bool, default=True) – In-place operations are not yet supported. Please set copy=True.
whiten (bool, default=False) – Not yet supported.
batch_size (int, optional) – Currently not needed and only added for API consistency and possible future extensions.
- Variables:
components (DNDarray of shape (n_components, n_features)) – Principal axes in feature space, representing the directions of maximum variance in the data. The components are sorted by `explained_variance_.
singular_values (DNDarray of shape (n_components,)) – The singular values corresponding to each of the selected components.
mean (DNDarray of shape (n_features,)) – Per-feature empirical mean, estimated from the training set.
n_components (int) – The estimated number of components.
n_samples_seen (int) – Number of samples processed so far.
- whiten = False
- n_components = None
- batch_size = None
- components_ = None
- singular_values_ = None
- mean_ = None
- n_components_ = None
- batch_size_ = None
- n_samples_seen_ = 0
- fit(path: str, chunk_size: int, dataset: str = 'DATA') Self[source]
Fit the IncrementalPCA model using data loaded in chunks from a HDF5 file.
This method processes data incrementally, loading chunks of data from a file and updating the PCA model iteratively. It is particularly useful for large datasets that cannot fit into memory.
- Parameters:
path (str) – Path to the file containing the dataset. The file must be in HDF5 format.
chunk_size (int) – Number of rows to load and process in each chunk. Must be smaller than or equal to the total number of rows in the dataset.
dataset (str, default="DATA") – Name of the dataset within the file to load.
- Returns:
The fitted IncrementalPCA instance.
- Return type:
Self
- Raises:
ValueError – If the file format is not HDF5. If chunk_size is larger than the number of rows in the dataset. If the number of columns is smaller than the number of processes.
- partial_fit(X: heat.DNDarray, y=None)[source]
One single step of incrementally building up the PCA. Input X is the current batch of data that needs to be added to the existing PCA.
- transform(X: heat.DNDarray) heat.DNDarray[source]
Apply dimensionality based on PCA to X.
- Parameters:
X (DNDarray of shape (n_samples, n_features)) – Data set to be transformed.
- inverse_transform(X: heat.DNDarray) heat.DNDarray[source]
Transform data back to its original space.
- Parameters:
X (DNDarray of shape (n_samples, n_components)) – Data set to be transformed back.
- _torch_matrix_diag(diagonal)
- class DMD(svd_solver: str | None = 'full', svd_rank: int | None = None, svd_tol: float | None = None)[source]
Bases:
heat.RegressionMixin,heat.BaseEstimatorDynamic Mode Decomposition (DMD), plain vanilla version with SVD-based implementation.
The time series of which DMD shall be computed must be provided as a 2-D DNDarray of shape (n_features, n_timesteps). Please, note that this deviates from Heat’s convention that data sets are handeled as 2-D arrays with the feature axis being the second axis.
- Parameters:
svd_solver (str, optional) – Specifies the algorithm to use for the singular value decomposition (SVD). Options are ‘full’ (default), ‘hierarchical’, and ‘randomized’.
svd_rank (int, optional) – The rank to which SVD shall be truncated. For ‘full’ SVD, svd_rank = None together with svd_tol = None (default) will result in no truncation. For svd_solver=’full’, at most one of svd_rank or svd_tol may be specified. For svd_solver=’hierarchical’, either svd_rank (rank to truncate to) or svd_tol (tolerance to truncate to) must be specified. For svd_solver=’randomized’, svd_rank must be specified and determines the the rank to truncate to.
svd_tol (float, optional) – The tolerance to which SVD shall be truncated. For ‘full’ SVD, svd_tol = None together with svd_rank = None (default) will result in no truncation. For svd_solver=’hierarchical’, either svd_tol (accuracy to truncate to) or svd_rank (rank to truncate to) must be specified. For svd_solver=’randomized’, svd_tol is meaningless and must be None.
- Variables:
svd_solver (str) – The algorithm used for the singular value decomposition (SVD).
svd_rank (int) – The rank to which SVD shall be truncated.
svd_tol (float) – The tolerance to which SVD shall be truncated.
rom_basis (DNDarray) – The reduced order model basis.
rom_transfer_matrix (DNDarray) – The reduced order model transfer matrix.
rom_eigenvalues (DNDarray) – The reduced order model eigenvalues.
rom_eigenmodes (DNDarray) – The reduced order model eigenmodes (“DMD modes”)
Notes
We follow the “exact DMD” method as described in [1], Sect. 2.2.
Please note that “rank” in the context of SVD always refers to the number of singular values/vectors to compute (i.e., “rank” refers to the mathematical rank of a matrix). This is completely different from the notion of “(MPI-)rank”, i.e., the ID given to a process, in a parallel MPI-application.
References
[1] J. L. Proctor, S. L. Brunton, and J. N. Kutz, “Dynamic Mode Decomposition with Control,” SIAM Journal on Applied Dynamical Systems, vol. 15, no. 1, pp. 142-161, 2016.
- svd_solver = 'full'
- svd_rank = None
- svd_tol = None
- rom_basis_ = None
- rom_transfer_matrix_ = None
- rom_eigenvalues_ = None
- rom_eigenmodes_ = None
- dmdmodes_ = None
- n_modes_ = None
- fit(X: heat.DNDarray) Self[source]
Fits the DMD model to the given data.
- Parameters:
X (DNDarray) – The time series data to fit the DMD model to. Must be of shape (n_features, n_timesteps).
- predict_next(X: heat.DNDarray, n_steps: int = 1) heat.DNDarray[source]
Predicts and returns the state(s) after n_steps-many time steps for given a current state(s).
- Parameters:
X (DNDarray) – The current state(s) for the prediction. Must have the same number of features as the training data, but can be batched for multiple current states, i.e., X can be of shape (n_features,) or (n_features, n_current_states). The output will have the same shape as the input.
n_steps (int, optional) – The number of steps to predict into the future. Default is 1, i.e., the next time step is predicted.
- predict(X: heat.DNDarray, steps: int | List[int]) heat.DNDarray[source]
Predics and returns future states given a current state(s) and returns them all as an array of size (n_steps, n_features).
This function avoids a time-stepping loop (i.e., repeated calls to ‘predict_next’) and computes the future states in one go. To do so, the number of future times to predict must be of moderate size as an array of shape (n_steps, self.n_modes_, self.n_modes_) must fit into memory. Moreover, it must be ensured that:
the array of initial states is not split or split along the batch axis (axis 1) and the feature axis is small (i.e., self.rom_basis_ is not split)
- Parameters:
X (DNDarray) – The current state(s) for the prediction. Must have the same number of features as the training data, but can be batched for multiple current states, i.e., X can be of shape (n_features,) or (n_current_states, n_features).
steps (int or List[int]) – if int: predictions at time step 0, 1, …, steps-1 are computed if List[int]: predictions at time steps given in the list are computed
- class DMDc(svd_solver: str | None = 'full', svd_rank: int | None = None, svd_tol: float | None = None)[source]
Bases:
heat.RegressionMixin,heat.BaseEstimatorDynamic Mode Decomposition with Control (DMDc), plain vanilla version with SVD-based implementation.
The time series of states and controls must be provided as 2-D DNDarrays of shapes (n_state_features, n_timesteps) and (n_control_features, n_timesteps), respectively. Please, note that this deviates from Heat’s convention that data sets are handeled as 2-D arrays with the feature axis being the second axis.
- Parameters:
svd_solver (str, optional) – Specifies the algorithm to use for the singular value decomposition (SVD). Options are ‘full’ (default), ‘hierarchical’, and ‘randomized’.
svd_rank (int, optional) – The rank to which SVD of the states shall be truncated. For ‘full’ SVD, svd_rank = None together with svd_tol = None (default) will result in no truncation. For svd_solver=’full’, at most one of svd_rank or svd_tol may be specified. For svd_solver=’hierarchical’, either svd_rank (rank to truncate to) or svd_tol (tolerance to truncate to) must be specified. For svd_solver=’randomized’, svd_rank must be specified and determines the the rank to truncate to.
svd_tol (float, optional) – The tolerance to which SVD of the states shall be truncated. For ‘full’ SVD, svd_tol = None together with svd_rank = None (default) will result in no truncation. For svd_solver=’hierarchical’, either svd_tol (accuracy to truncate to) or svd_rank (rank to truncate to) must be specified. For svd_solver=’randomized’, svd_tol is meaningless and must be None.
- Variables:
svd_solver (str) – The algorithm used for the singular value decomposition (SVD).
svd_rank (int) – The rank to which SVD shall be truncated.
svd_tol (float) – The tolerance to which SVD shall be truncated.
rom_basis (DNDarray) – The reduced order model basis.
rom_transfer_matrix (DNDarray) – The reduced order model transfer matrix.
rom_control_matrix (DNDarray) – The reduced order model control matrix.
rom_eigenvalues (DNDarray) – The reduced order model eigenvalues.
rom_eigenmodes (DNDarray) – The reduced order model eigenmodes (“DMD modes”)
Notes
We follow the approach described in [1], Sects. 3.3 and 3.4. In the case that svd_rank is prescribed, the rank of the SVD of the full system matrix is set to svd_rank + n_control_features; cf. https://github.com/dynamicslab/pykoopman for the same approach.
Please note that “rank” in the context of SVD always refers to the number of singular values/vectors to compute (i.e., “rank” refers to the mathematical rank of a matrix). This is completely different from the notion of “(MPI-)rank”, i.e., the ID given to a process, in a parallel MPI-application.
References
[1] J. L. Proctor, S. L. Brunton, and J. N. Kutz, “Dynamic Mode Decomposition with Control,” SIAM Journal on Applied Dynamical Systems, vol. 15, no. 1, pp. 142-161, 2016.
- svd_solver = 'full'
- svd_rank = None
- svd_tol = None
- rom_basis_ = None
- rom_transfer_matrix_ = None
- rom_control_matrix_ = None
- rom_eigenvalues_ = None
- rom_eigenmodes_ = None
- dmdmodes_ = None
- n_modes_ = None
- n_modes_system_ = None
- fit(X: heat.DNDarray, C: heat.DNDarray) Self[source]
Fits the DMD model to the given data.
- predict(X: heat.DNDarray, C: heat.DNDarray) heat.DNDarray[source]
Predicts and returns future states given the current state(s)
Xand control trajectoryC.- Parameters:
X (DNDarray) – The current state(s) for the prediction. Must have the same number of features as the training data, but can be batched for multiple current states, i.e., X can be of shape (n_state_features,) or (n_batch, n_state_features).
C (DNDarray) – The control trajectory for the prediction. Must have the same number of control features as the training data, i.e., C must be of shape (n_control_features,) –for a single time step– or (n_control_features, n_timesteps).