POD

Module for Proper Orthogonal Decomposition (POD). Three different methods can be employed: Truncated Singular Value Decomposition, Truncated Randomized Singular Value Decomposition, Truncated Singular Value Decomposition via correlation matrix.

class POD(method='svd', rank=-1, subspace_iteration=1, omega_rank=0, save_memory=False)[source]

Bases: Reduction

Perform the Proper Orthogonal Decomposition.

Parameters:
  • method ({'svd', 'randomized_svd', 'correlation_matrix'}) – the implementation to use for the computation of the POD modes. Default is ‘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. The rank parameter is available using all the available methods.

  • subspace_iteration (int) – the number of subspace iteration in the randomized svd. It is available only using the ‘randomized_svd’ method. Default value is 1.

  • omega_rank (int) – the number of columns of the Omega random matrix. If set to 0, the number of columns is equal to twice the rank (if it has explicitly passed as integer) or twice the number of input snapshots. Default value is 0. It is available only using the ‘randomized_svd’ method.

  • save_memory (bool) – reduce the usage of the memory, despite an higher number of operations. It is available only using the ‘correlation_matrix’ method. Default value is False.

Example:
>>> pod = POD().fit(snapshots)
>>> reduced_snapshots = pod.reduce(snapshots)
>>> # Other possible constructors are ...
>>> pod = POD('svd')
>>> pod = POD('svd', rank=20)
>>> pod = POD('randomized_svd', rank=-1)
>>> pod = POD('randomized_svd', rank=0, subspace_iteration=3,
              omega_rank=10)
>>> pod = POD('correlation_matrix', rank=10, save_memory=False)
_abc_impl = <_abc._abc_data object>
_corrm(X)[source]

Truncated POD calculated with correlation matrix.

Parameters:

X (numpy.ndarray) – the matrix to decompose.

Returns:

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

Return type:

numpy.ndarray, numpy.ndarray, numpy.ndarray

_rsvd(X)[source]

Truncated randomized Singular Value Decomposition.

Computes an approximate SVD using randomized algorithms for efficiency.

Parameters:

X (numpy.ndarray) – The matrix to decompose.

Returns:

Tuple of (truncated left-singular vectors, truncated singular values).

Return type:

tuple(numpy.ndarray, numpy.ndarray)

References:

Finding structure with randomness: probabilistic algorithms for constructing approximate matrix decompositions. N. Halko, P. G. Martinsson, J. A. Tropp.

_svd(X)[source]

Truncated Singular Value Decomposition.

Parameters:

X (numpy.ndarray) – the matrix to decompose.

Returns:

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

Return type:

numpy.ndarray, numpy.ndarray, numpy.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.

_truncation(X, s)[source]

Return the number of modes to select according to the rank parameter. See POD.__init__ for further info.

Parameters:
Returns:

the number of modes

Return type:

int

expand(X)[source]

Projects a reduced to full order solution.

Type:

numpy.ndarray

Note

Same as inverse_transform. Kept for backward compatibility.

fit(X)[source]

Create the reduced space for the given snapshots using POD.

Computes the POD modes and singular values using the specified method.

Parameters:

X (numpy.ndarray) – The input snapshots matrix (stored by column).

Returns:

self

inverse_transform(X)[source]

Projects a reduced to full order solution.

Type:

numpy.ndarray

property modes

The POD modes.

Type:

numpy.ndarray

reduce(X)[source]

Reduces the given snapshots.

Parameters:

X (numpy.ndarray) – the input snapshots matrix (stored by column).

Note

Same as transform. Kept for backward compatibility.

property singular_values

The singular values

Type:

numpy.ndarray

transform(X)[source]

Reduces the given snapshots.

Parameters:

X (numpy.ndarray) – the input snapshots matrix (stored by column).