Compressed DMD

Derived module from dmdbase.py for compressed dmd. As a reference consult this work by Erichson, Brunton and Kutz: https://doi.org/10.1007/s11554-016-0655-2

class CDMD(svd_rank=0, tlsq_rank=0, compression_matrix='uniform', opt=False, rescale_mode=None, forward_backward=False, sorted_eigs=False, tikhonov_regularization=None)[source]

Bases: DMDBase

Compressed Dynamic Mode Decomposition.

Compute the dynamic mode decomposition after the multiplication of the snapshots matrix by the compression_matrix, in order to compress the input data. It is possible use a custom matrix for the compression or chose between the preconstructed matrices. Available values for compression_matrix are:

  • ‘normal’: the matrix C with dimension (nsnaps, ndim) is randomly generated with normal distribution with mean equal to 0.0 and standard deviation equal to 1.0;

  • ‘uniform’: the matrix C with dimension (nsnaps, ndim) is randomly generated with uniform distribution between 0 and 1;

  • ‘sparse’: the matrix C with dimension (nsnaps, ndim) is a random sparse matrix;

  • ‘sample’: the matrix C with dimension (nsnaps, ndim) where each row contains an element equal to 1 and all the other elements are null.

Parameters:
  • svd_rank (int or float) – the rank for the truncation; If 0, the method computes the optimal rank and uses it for truncation; if positive interger, the method uses the argument for the truncation; if float between 0 and 1, the rank is the number of the biggest singular values that are needed to reach the ‘energy’ specified by svd_rank; if -1, the method does not compute truncation.

  • tlsq_rank (int) – rank truncation computing Total Least Square. Default is 0, that means TLSQ is not applied.

  • compression_matrix ({'linear', 'sparse', 'uniform', 'sample'} or numpy.ndarray) – the matrix that pre-multiplies the snapshots matrix in order to compress it; if compression_matrix is a numpy.ndarray, its dimension must be (nsnaps, ndim). Default value is ‘uniform’.

  • opt (bool or int) – argument to control the computation of DMD modes amplitudes. See DMDBase. Default is False.

  • rescale_mode ({'auto'} or None or numpy.ndarray) – Scale Atilde as shown in 10.1016/j.jneumeth.2015.10.010 (section 2.4) before computing its eigendecomposition. None means no rescaling, ‘auto’ means automatic rescaling using singular values, otherwise the scaling factors.

  • forward_backward (bool) – If True, the low-rank operator is computed like in fbDMD (reference: https://arxiv.org/abs/1507.02264). Default is False.

  • sorted_eigs ({'real', 'abs'} or False) – Sort eigenvalues (and modes/dynamics accordingly) by magnitude if sorted_eigs=’abs’, by real part (and then by imaginary part to break ties) if sorted_eigs=’real’.

  • sorted_eigs – Sort eigenvalues (and modes/dynamics accordingly) by magnitude if sorted_eigs=’abs’, by real part (and then by imaginary part to break ties) if sorted_eigs=’real’. Default: False.

  • tikhonov_regularization (int or float) – Tikhonov parameter for the regularization. If None, no regularization is applied, if float, it is used as the \lambda tikhonov parameter.

_compress_snapshots()[source]

Private method that compresses the snapshots matrix by pre-multiplying it by the chosen compression_matrix.

Returns:

the compressed snapshots.

Return type:

numpy.ndarray

property compression_matrix

The compression matrix

fit(X)[source]

Compute the Dynamic Modes Decomposition to the input data.

Parameters:

X (numpy.ndarray or iterable) – the input snapshots.