Sparsity-promoting DMD¶
Derived module from dmdbase.py for sparsity-promoting DMD.
- class SpDMD(svd_rank=0, tlsq_rank=0, exact=True, opt=False, rescale_mode=None, forward_backward=False, sorted_eigs=False, abs_tolerance=1e-06, rel_tolerance=0.0001, max_iterations=10000, rho=1, gamma=10, verbose=True, enforce_zero=True, release_memory=True, zero_absolute_tolerance=1e-12)[source]
Bases:
DMD
Sparsity-Promoting Dynamic Mode Decomposition. Promotes solutions having an high number of amplitudes set to zero (i.e. sparse solutions). Reference: 10.1063/1.4863670
- 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.
exact (bool) – flag to compute either exact DMD or projected DMD. Default is True.
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’. Default: False.
abs_tolerance (float) – Controls the convergence of ADMM. See
_loop_condition()
for more details.rel_tolerance (float) – Controls the convergence of ADMM. See
_loop_condition()
for more details.max_iterations (int) – The maximum number of iterations performed by ADMM, after that the algorithm is stopped.
rho (float) – Controls the convergence of ADMM. For a reference on the optimal value for rho see 10.1109/TAC.2014.2354892 or 10.3182/20120914-2-US-4030.00038.
gamma (float) – Controls the level of “promotion” assigned to sparse solution. Increasing gamma will result in an higher number of zero-amplitudes.
verbose (bool) – If False, the information provided by SpDMD (like the number of iterations performed by ADMM) are not shown.
enforce_zero (bool) – If True the DMD amplitudes which should be set to zero according to the solution of ADMM are manually set to 0 (since we solve a sparse linear system to find the optimal vector of DMD amplitudes very small terms may survive in some cases).
release_memory – If True the intermediate matrices computed by the algorithm are deleted after the termination of a call to
fit()
.
- _find_zero_amplitudes()[source]
Use ADMM to find which amplitudes (i.e. their position in the DMD amplitudes array) which can be set to zero according to the given parameters. Note that this method does not compute amplitudes, but only which amplitudes are to be set to 0. Optimal amplitudes should be computed separately afterwards (see
_find_sparsity_promoting_amplitudes()
). :return np.ndarray: A boolean vector whose True items correspond toamplitudes which should be set to 0.
- _loop_condition(alpha, beta, lmbd, old_beta)[source]
Check whether ADMM can stop now, or should perform another iteration. :param np.ndarray alpha: Current value of \alpha_k (vector
of DMD amplitudes).
- Parameters:
beta (np.ndarray) – Current value of \beta_k (vector of non-zero amplitudes).
lmbd (np.ndarray) – Current value of \lambda_k (vector of Lagrange multipliers).
old_beta (np.ndarray) – Old value of \beta_{k-1} (vector of non-zero amplitudes).
- Return bool:
True if ADMM can stop now, False otherwise.
- _optimal_amplitudes(zero_amplitudes)[source]
Find the optimal DMD amplitudes with the constraint that the given indexes should be set to 0. :param np.ndarray zero_amplitudes: Boolean vector. :return np.ndarray: Vector of optimal DMD amplitudes. Amplitudes at
indexes corresponding to True indexes in zero_amplitudes are set to 0.
- _update(beta, lmbd)[source]
Operate an entire step of ADMM. :param np.ndarray beta: Current value of \beta_k (vector of
non-zero amplitudes).
- Parameters:
lmbd (np.ndarray) – Current value of \lambda_k (vector of Lagrande multipliers).
- Returns:
A tuple containing the updated values \alpha_{k+1},\beta_{k+1},\lambda_{k+1} (in this order).
- Return type:
tuple
- _update_alpha(beta, lmbd)[source]
Update the vector \alpha_k of DMD amplitudes. :param np.ndarray beta: Current value of \beta_k (vector of
non-zero amplitudes).
- Parameters:
lmbd (np.ndarray) – Current value of \lambda_k (vector of Lagrande multipliers).
- Returns:
The updated value \alpha_{k+1}.
- Return type:
np.ndarray
- _update_beta(alpha, lmbd)[source]
Update the vector \beta of non-zero amplitudes. :param np.ndarray alpha: Updated value of \alpha_{k+1} (vector
of DMD amplitudes).
- Parameters:
lmbd (np.ndarray) – Current value of \lambda_k (vector of Lagrange multipliers).
- Returns:
The updated value \beta_{k+1}.
- Return type:
np.ndarray
- _update_lagrangian(alpha, beta, lmbd)[source]
Update the vector \lambda of Lagrange multipliers. :param np.ndarray alpha: Updated value of \alpha_{k+1} (vector
of DMD amplitudes).
- Parameters:
beta (np.ndarray) – Updated value of \beta_{k+1} (vector of non-zero amplitudes).
lmbd (np.ndarray) – Current value of \lambda_k (vector of Lagrange multipliers).
- Returns:
The updated value \lambda_{k+1}.
- Return type:
np.ndarray
- fit(X)[source]
Compute the Dynamic Modes Decomposition of the input data. :param X: the input snapshots. :type X: numpy.ndarray or iterable