DMD Operator

class DMDOperator(svd_rank, exact, forward_backward, rescale_mode, sorted_eigs, tikhonov_regularization)[source]

Bases: object

Dynamic Mode Decomposition standard operator class. Non-standard ways of computing the low-rank Atilde operator should be coded into subclasses of this class.

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.

  • exact (bool) – flag to compute either exact DMD or projected DMD. 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’. 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.

property Lambda
_compute_eigenquantities()[source]

Private method that computes eigenvalues and eigenvectors of the low-dimensional operator, scaled according to self._rescale_mode.

_compute_modes(Y, U, Sigma, V)[source]

Private method that computes eigenvalues and eigenvectors of the high-dimensional operator (stored in self.modes and self.Lambda).

Parameters
_least_square_operator(U, s, V, Y)[source]

Private method that computes the lowrank operator from the singular value decomposition of matrix X and the matrix Y.

\mathbf{\tilde{A}} = \mathbf{U}^* \mathbf{Y} \mathbf{X}^\dagger \mathbf{U} = \mathbf{U}^* \mathbf{Y} \mathbf{V} \mathbf{S}^{-1}

Parameters
  • U (numpy.ndarray) – 2D matrix that contains the left-singular vectors of X, stored by column.

  • s (numpy.ndarray) – 1D array that contains the singular values of X.

  • V (numpy.ndarray) – 2D matrix that contains the right-singular vectors of X, stored by row.

  • Y (numpy.ndarray) – input matrix Y.

Returns

the lowrank operator

Return type

numpy.ndarray

property as_numpy_array
compute_operator(X, Y)[source]

Compute the low-rank operator.

Parameters
  • X (numpy.ndarray) – matrix containing the snapshots x0,..x{n-1} by column.

  • Y (numpy.ndarray) – matrix containing the snapshots x1,..x{n} by column.

Returns

the (truncated) left-singular vectors matrix, the (truncated) singular values array, the (truncated) right-singular vectors matrix of X.

Return type

numpy.ndarray, numpy.ndarray, numpy.ndarray

property eigenvalues
property eigenvectors
property modes
property shape

Shape of the operator