Kernelized extended DMD

Derived module from dmdbase.py for kernelized extended dmd.

References: - M. O. Williams, C. W. Rowley, and I. G. Kevrekidis, A kernel-based method for data-driven koopman spectral analysis, J. Comput. Dynam., 2 (2015), pp. 247-265 - M. O. Williams, I. G. Kevrekidis, and C. W. Rowley, A data-driven approximation of the koopman operator: extending dynamic mode decomposition, J. Nonlinear Sci., 25 (2015), pp. 1307-1346

class EDMD(svd_rank=-1, tlsq_rank=0, opt=False, kernel_metric='linear', kernel_params=None)[source]

Bases: pydmd.dmd.DMD

Kernelized Extended DMD.

Parameters
  • svd_rank (int or float) – the rank for the truncation; If positive integer, 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, meaning no truncation.

  • opt (bool or int) – If True, amplitudes are computed like in optimized DMD (see _compute_amplitudes() for reference). If False, amplitudes are computed following the standard algorithm. If opt is an integer, it is used as the (temporal) index of the snapshot used to compute DMD modes amplitudes (following the standard algorithm). The reconstruction will generally be better in time instants near the chosen snapshot; however increasing opt may lead to wrong results when the system presents small eigenvalues. For this reason a manual selection of the number of eigenvalues considered for the analyisis may be needed (check svd_rank). Also setting svd_rank to a value between 0 and 1 may give better results. Default is False.

  • kernel_metric (str) – the kernel function to apply. Supported kernel metrics include “additive_chi2”, “chi2”, “linear”, “poly”, “polynomial”, “rbf”, “laplacian”, “sigmoid”, and “cosine”.

  • kernel_params (dict) – additional parameters for the sklearn.metrics.pairwise_kernels function, including kernel-specific function parameters.

static _test_kernel_inputs(kernel_metric, kernel_params)[source]

Helper function that uses a dummy array of data in order to call sklearn.metrics.pairwise_kernels using the user-given kernel parameters. Ensures that the given kernel parameters produce valid kernel matrices.

eigenfunctions(x)[source]
Parameters

x (numpy.ndarray) – array from the original snapshot domain denoting the point at which to compute the EDMD eigenfunctions.

Returns

array containing all r EDMD eigenfunctions evaluated at the given domain point, where r denotes the rank of the EDMD fit.

Return type

numpy.ndarray

class EDMDOperator(svd_rank, kernel_metric, kernel_params)[source]

Bases: pydmd.dmdoperator.DMDOperator

DMD operator for kernelized extended DMD.

Parameters
  • svd_rank (int or float) – the rank for the truncation; If positive integer, 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.

  • kernel_metric (str) – the kernel function to apply. Supported kernel metrics include “additive_chi2”, “chi2”, “linear”, “poly”, “polynomial”, “rbf”, “laplacian”, “sigmoid”, and “cosine”.

  • kernel_params (dict) – additional parameters for the sklearn.metrics.pairwise_kernels function, including kernel-specific function parameters.

_compute_feature_matrix_svd(Gramian)[source]

Helper function that computes the eigendecomposition of the Gramian in order to obtain the singular vectors Q and singular values Sigma of the corresponding feature embedding matrix.

Parameters

Gramian (numpy.ndarray) – the Gramian matrix.

Returns

the (truncated) left-singular vectors matrix of the extended feature matrix and the (truncated) singular values matrix of the extended feature matrix.

Return type

numpy.ndarray, numpy.ndarray

_compute_rank_from_svals(s)[source]

Rank computation for the truncated Singular Value Decomposition given only singular values.

Parameters

s (numpy.ndarray) – the singular values.

Returns

the computed rank truncation.

Return type

int

property as_numpy_array

EDMD does not explicitly compute the forward operator or its reduction. Doing so poses risks of prohibitively-large computations…

compute_operator(X, Y)[source]

Compute and store the kernelized EDMD diagnostics.

Parameters
  • X (numpy.ndarray) – matrix containing the right-hand side snapshots by column.

  • Y (numpy.ndarray) – matrix containing the left-hand side snapshots by column.

Returns

the (truncated) left-singular vectors matrix of the extended feature matrix, the (truncated) singular values matrix of the extended feature matrix, and the eigenvectors of the tranformed forward operator K_hat.

Return type

numpy.ndarray, numpy.ndarray, numpy.ndarray

property svd_vals
Returns

the singular values of the extended feature matrix.

Return type

numpy.ndarray