Utilities

Utilities module.

compute_rank(X: numpy.ndarray, svd_rank: numbers.Number = 0) → int[source]

Rank computation for the truncated Singular Value Decomposition.

Parameters
  • X (np.ndarray) – the matrix to decompose.

  • 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. Default is 0.

Returns

the computed rank truncation.

Return type

int

References: Gavish, Matan, and David L. Donoho, The optimal hard threshold for singular values is, IEEE Transactions on Information Theory 60.8 (2014): 5040-5053.

compute_svd(X: numpy.ndarray, svd_rank: numbers.Number = 0) → pydmd.utils.SVD[source]

Truncated Singular Value Decomposition.

Parameters
  • X (np.ndarray) – the matrix to decompose.

  • 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. Default is 0.

Returns

the truncated left-singular vectors matrix, the truncated singular values array, the truncated right-singular vectors matrix.

Return type

NamedTuple(“SVD”, [(‘U’, np.ndarray), (‘s’, np.ndarray), (‘V’, np.ndarray)])

References: Gavish, Matan, and David L. Donoho, The optimal hard threshold for singular values is, IEEE Transactions on Information Theory 60.8 (2014): 5040-5053.

compute_tlsq(X: numpy.ndarray, Y: numpy.ndarray, tlsq_rank: int) → pydmd.utils.TLSQ[source]

Compute Total Least Square.

Parameters
  • X (np.ndarray) – the first matrix;

  • Y (np.ndarray) – the second matrix;

  • tlsq_rank (int) – the rank for the truncation; If 0, the method does not compute any noise reduction; if positive number, the method uses the argument for the SVD truncation used in the TLSQ method.

Returns

the denoised matrix X, the denoised matrix Y

Return type

NamedTuple(“TLSQ”, [(‘X_denoised’, np.ndarray), (‘Y_denoised’, np.ndarray)])

References: https://arxiv.org/pdf/1703.11004.pdf https://arxiv.org/pdf/1502.03854.pdf

pseudo_hankel_matrix(X: numpy.ndarray, d: int) → numpy.ndarray[source]

Arrange the snapshot in the matrix X into the (pseudo) Hankel matrix. The attribute d controls the number of snapshot from X in each snapshot of the Hankel matrix.

Example
>>> a = np.array([[1, 2, 3, 4, 5]])
>>> _hankel_pre_processing(a, d=2)
array([[1, 2, 3, 4],
       [2, 3, 4, 5]])
>>> _hankel_pre_processing(a, d=4)
array([[1, 2],
       [2, 3],
       [3, 4],
       [4, 5]])
>>> a = np.array([1,2,3,4,5,6]).reshape(2,3)
array([[1, 2, 3],
       [4, 5, 6]])
>>> _hankel_pre_processing(a, d=2)
array([[1, 2],
       [4, 5],
       [2, 3],
       [5, 6]])