VarProDMD: Variable Projection for DMD¶
Variable Projection for DMD. Reformulation of original paper (https://epubs.siam.org/doi/abs/10.1137/M1124176) s.t. sparse matrix computation is substiuted by outer products. Further the optimization is reformulated s.t. SciPy’s nonlinear least squares optimizer can handle “complex” parameters.
- Default optimizer arguments:
- OPT_DEF_ARGS =
{“method”: “trf”, “tr_solver”: “exact”, “loss”: “linear”, “x_scale”: “jac”, “gtol”: 1e-8, “xtol”: 1e-8, “ftol”: 1e-8}
- class VarProDMD(svd_rank: float | int = 0, exact: bool = False, sorted_eigs: bool | str = False, compression: float = 0.0, optargs: Dict[str, Any] | None = None)[source]
Bases:
DMDBase
Variable Projection for DMD using SciPy’s nonlinear least squares solver. The original problem is reformulated s.t. complex residuals and Jacobians, which are used by the Levenberg-Marquardt algorithm, are transormed to real numbers. Further simplifications (outer products) avoids using sparse matrices.
VarProDMD constructor.
- Parameters:
svd_rank (Union[float, int], optional) – Desired SVD rank. If rank r = 0, the optimal rank is determined automatically. If rank is a float s.t. 0 < r < 1, the cumulative energy of the singular values is used to determine the optimal rank. If rank is an integer and r > 0, the desired rank is used iff possible. Defaults to 0.
exact (bool, optional) – Perform variable projection in low dimensional space if exact=False. Else the optimization is performed in the original space. Defaults to False.
sorted_eigs (Union[bool, str], optional) – Sort eigenvalues. If sorted_eigs=True, the variance of the absolute values of the complex eigenvalues \left(\sqrt{\omega_i \cdot \bar{\omega}_i}\right), the variance absolute values of the real parts \left|\Re\{{\omega_i}\}\right| and the variance of the absolute values of the imaginary parts \left|\Im\{{\omega_i}\}\right| is computed. The eigenvalues are then sorted according to the highest variance (from highest to lowest). If sorted_eigs=False, no sorting is performed. If the parameter is a string and set to sorted_eigs=’auto’, the eigenvalues are sorted accoring to the variances of previous mentioned quantities. If sorted_eigs=’real’ the eigenvalues are sorted w.r.t. the absolute values of the real parts of the eigenvalues (from highest to lowest). If sorted_eigs=’imag’ the eigenvalues are sorted w.r.t. the absolute values of the imaginary parts of the eigenvalues (from highest to lowest). If sorted_eigs=’abs’ the eigenvalues are sorted w.r.t. the magnitude of the eigenvalues \left(\sqrt{\omega_i \cdot \bar{\omega}_i}\right) (from highest to lowest). Defaults to False.
compression (float, optional) – If libary compression c = 0, all samples are used. If 0 < c < 1, the best fitting \lfloor \left(1 - c\right)m\rfloor samples are selected.
optargs (Dict[str, Any], optional) – Arguments for ‘least_squares’ optimizer. If set to None, OPT_DEF_ARGS are used as default parameters. Defaults to None.
- property dynamics
Get the time evolution of each mode.
- Returns:
matrix that contains all the time evolution, stored by row.
- Return type:
- fit(X: ndarray, time: ndarray) object [source]
Fit the eigenvalues, modes and eigenfunctions/amplitudes to measurements.
- Parameters:
X (np.ndarray) – Measurements \boldsymbol{X} \in \mathbb{C}^{n \times m}.
time (np.ndarray) – 1d array of timestamps where measurements were taken.
- Returns:
VarProDMD instance.
- Return type:
- forecast(time: ndarray) ndarray [source]
Forecast measurements at given timestamps time.
- Parameters:
time (np.ndarray) – Desired times for forcasting as 1d array.
- Raises:
ValueError – If method fit(X, time) was not called.
- Returns:
Predicted measurements \hat{\boldsymbol{X}}.
- Return type:
np.ndarray
- property frequency
Get the amplitude spectrum.
- Returns:
the array that contains the frequencies of the eigenvalues.
- Return type:
- property growth_rate
Get the growth rate values relative to the modes.
- Returns:
the Floquet values
- Return type:
- property opt_stats: OptimizeResult
Return optimization statistics of the Variable Projection optimization.
- Raises:
ValueError – ValueError is raised if method fit(X, time) was not called.
- Returns:
Optimization results including optimal weights (continuous eigenvalues) and number of iterations.
- Return type:
OptimizeResult
- property selected_samples: ndarray
Return indices for creating the library.
- Raises:
ValueError – ValueError is raised if method fit(X, time) was not called.
- Returns:
Indices of the selected samples. If no compression was performed \left(c = 0\right), all indices are returned, else indices of the library selection scheme using QR-Decomposition with column pivoting.
- Return type:
np.ndarray
- property ssr: float
Compute the square root of sum squared residual (SSR) taken from https://link.springer.com/article/10.1007/s10589-012-9492-9. The SSR gives insight w.r.t. signal qualities. A low SSR is desired. If SSR is high the model may be inaccurate.
- Raises:
ValueError – ValueError is raised if method fit(X, time) was not called.
- Returns:
SSR.
- Return type:
- class VarProOperator(svd_rank: float | int, exact: bool, sorted_eigs: bool | str, compression: float, optargs: Dict[str, Any])[source]
Bases:
DMDOperator
Variable Projection Operator for DMD.
VarProOperator constructor.
- Parameters:
svd_rank (Union[float, int]) – Desired SVD rank. If rank r = 0, the optimal rank is determined automatically. If rank is a float s.t. 0 < r < 1, the cumulative energy of the singular values is used to determine the optimal rank. If rank is an integer and r > 0, the desired rank is used iff possible.
exact (bool) – Perform computations in original state space if exact=True, else perform optimization in projected (low dimensional) space.
sorted_eigs (Union[bool, str]) – Sort eigenvalues. If sorted_eigs=True, the variance of the absolute values of the complex eigenvalues \sqrt{\omega_i \cdot \bar{\omega}_i}, the variance absolute values of the real parts \left|\Re\{{\omega_i}\}\right| and the variance of the absolute values of the imaginary parts \left|\Im\{{\omega_i}\}\right| is computed. The eigenvalues are then sorted according to the highest variance (from highest to lowest). If sorted_eigs=False, no sorting is performed. If the parameter is a string and set to sorted_eigs=’auto’, the eigenvalues are sorted accoring to the variances of previous mentioned quantities. If sorted_eigs=’real’ the eigenvalues are sorted w.r.t. the absolute values of the real parts of the eigenvalues (from highest to lowest). If sorted_eigs=’imag’ the eigenvalues are sorted w.r.t. the absolute values of the imaginary parts of the eigenvalues (from highest to lowest). If sorted_eigs=’abs’ the eigenvalues are sorted w.r.t. the magnitude of the eigenvalues \left(\sqrt{\omega_i \cdot \bar{\omega}_i}\right) (from highest to lowest).
compression (float) – If libary compression c = 0, all samples are used. If 0 < c < 1, the best fitting \lfloor \left(1 - c\right)m\rfloor samples are selected.
optargs (Dict[str, Any]) – Arguments for ‘least_squares’ optimizer. Use OPT_DEF_ARGS as starting point.
- compute_operator(X: ndarray, Y: ndarray) Tuple[ndarray, OptimizeResult, ndarray] [source]
Perform Variable Projection for DMD using SciPy’s nonlinear least squares optimizer.
- Parameters:
X (np.ndarray) – Measurement \boldsymbol{X} \in \mathbb{C}^{n \times m}
Y (np.ndarray) – 1d array of timestamps where individual measurements \boldsymbol{x}_i \in \mathbb{C}^n where taken.
- Raises:
ValueError – If sorted_eigs from constructor was set to an invalid string.
- Returns:
Eigenfunctions/amplitudes \boldsymbol{\varphi}^{m-1}, OptimizeResult from SciPy’s optimizer (optimal parameters and statistics), indices of selected measurements. If no compression (c = 0) is used, all indices are returned, else the indices of the selected samples are used.
- Return type:
Tuple[np.ndarray, OptimizeResult, np.ndarray]