Hankel DMD

Derived module from dmdbase.py for hankel dmd.

Reference: - H. Arbabi, I. Mezic, Ergodic theory, dynamic mode decomposition, and computation of spectral properties of the Koopman operator. SIAM Journal on Applied Dynamical Systems, 2017, 16.4: 2096-2126.

class HankelDMD(svd_rank=0, tlsq_rank=0, exact=False, opt=False, rescale_mode=None, forward_backward=False, d=1, sorted_eigs=False, reconstruction_method='first', tikhonov_regularization=None)[source]

Bases: pydmd.dmdbase.DMDBase

Hankel Dynamic Mode Decomposition

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 no truncation.

  • exact (bool) – flag to compute either exact DMD or projected DMD. Default is False.

  • 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.

  • d (int) – the new order for spatial dimension of the input snapshots. Default is 1.

  • 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.

  • reconstruction_method ({'first', 'mean'} or array-like) – Method used to reconstruct the snapshots of the dynamical system from the multiple versions available due to how HankelDMD is conceived. If ‘first’ (default) the first version available is selected (i.e. the nearest to the 0-th row in the augmented matrix). If ‘mean’ we compute the element-wise mean. If reconstruction_method is an array of float values we compute the weighted average (for each snapshots) using the given values as weights (the number of weights must be equal to d).

_hankel_first_occurrence(time)[source]

For a given t such that there is k \in \mathbb{N} such that t = t_0 + k dt, return the index of the first column in Hankel pseudo matrix (see also _pseudo_hankel_matrix()) which contains the snapshot corresponding to t.

Parameters

time – The time corresponding to the requested snapshot.

Returns

The index of the first appeareance of time in the columns of Hankel pseudo matrix.

Return type

int

_update_sub_dmd_time()[source]

Update the time dictionaries (dmd_time and original_time) of the auxiliary DMD instance HankelDMD._sub_dmd after an update of the time dictionaries of the time dictionaries of this instance of the higher level instance of HankelDMD.

property amplitudes

Get the coefficients that minimize the error between the original system and the reconstructed one. For futher information, see dmdbase._compute_amplitudes.

Returns

the array that contains the amplitudes coefficient.

Return type

numpy.ndarray

property eigs

Get the eigenvalues of A tilde.

Returns

the eigenvalues from the eigendecomposition of atilde.

Return type

numpy.ndarray

fit(X)[source]

Compute the Dynamic Modes Decomposition to the input data.

Parameters

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

property ho_snapshots

Get the time-delay data matrix.

Returns

the matrix that contains the time-delayed data.

Return type

numpy.ndarray

property modes

Get the matrix containing the DMD modes, stored by column.

Returns

the matrix containing the DMD modes.

Return type

numpy.ndarray

property modes_activation_bitmask

Get the bitmask which controls which DMD modes are enabled at the moment in this DMD instance.

The DMD instance must be fitted before this property becomes valid. After fit() is called, the defalt value of modes_activation_bitmask is an array of True values of the same shape of amplitudes().

The array returned is read-only (this allow us to react appropriately to changes in the bitmask). In order to modify the bitmask you need to set the field to a brand-new value (see example below).

Example:

>>> # this is an error
>>> dmd.modes_activation_bitmask[[1,2]] = False
ValueError: assignment destination is read-only
>>> tmp = np.array(dmd.modes_activation_bitmask)
>>> tmp[[1,2]] = False
>>> dmd.modes_activation_bitmask = tmp
Returns

The DMD modes activation bitmask.

Return type

numpy.ndarray

property operator

Get the instance of DMDOperator.

Returns

the instance of DMDOperator

Return type

DMDOperator

property reconstructed_data

Get the reconstructed data.

Returns

the matrix that contains the reconstructed snapshots.

Return type

numpy.ndarray

property svd_rank