PODBasis#

class PODBasis(num_vectors: int = None, svdval_threshold: float = None, cumulative_energy: float = None, residual_energy: float = None, projection_error: float = None, max_vectors: int = None, svdsolver: str = 'dense', weights: ndarray = None, name: str = None, **svdsolver_options)[source]#

Proper othogonal decomposition basis, consisting of the principal left singular vectors of a collection of states.

\[\text{svd}(\Q) = \bfPhi\bfSigma\bfPsi\trp \qquad\Longrightarrow\qquad \text{pod}(\Q, r) = \bfPhi_{:,:r}\]

The low-dimensional approximation is linear, see LinearBasis. Here, \(\Q\in\mathbb{R}^{n\times k}\) is a collection of states, the only argument of fit().

The POD basis entries \(\Vr = \bfPhi_{:,:r}\in\RR^{n\times r}\) are always orthogonal, i.e., \(\Vr\trp\Vr = \I\). If a weight matrix \(\W\) is specified, a weighted SVD is computed so that \(\Vr\trp\W\Vr = \I\).

The number of left singular vectors \(r\) is the dimension of the reduced state and is set by specifying exactly one of the constructor arguments num_vectors, svdval_threshold, residual_energy, cumulative_energy, or projection_error. Once the basis entries are set by calling fit(), the reduced state dimension \(r\) can be updated by calling set_dimension().

The POD singular values, which are used to select \(r\), are the diagonals of \(\bfSigma\) and are denoted \(\sigma_1,\ldots,\sigma_k\).

Parameters
num_vectorsint

Set the reduced state dimension \(r\) to num_vectors.

svdval_thresholdfloat

Choose \(r\) as the number of normalized POD singular values that are greater than the given threshold, i.e., \(\sigma_{i}/\sigma_{1} \ge\) svdval_threshold for \(i=1,\ldots,r\).

cumulative_energyfloat

Choose \(r\) as the smallest integer such that \(\sum_{i=1}^{r}\sigma_i^2\big/\sum_{j=1}^{k}\sigma_j^2\) is greater than or equal to cumulative_energy.

residual_energyfloat

Choose \(r\) as the smallest integer such that \(\sum_{i=r+1}^k\sigma_i^2\big/\sum_{j=1}^k\sigma_j^2\) is less than or equal to residual_energy.

projection_errorfloat

Choose \(r\) as the smallest integer such that \(\|\Q - \Vr\Vr\trp\Q\|_F \big/ \|\Q\|_F\) is less than or equal to projection_error.

max_vectorsint

Maximum number of POD basis vectors to store. After calling fit(), the reduced_state_dimension can be increased up to max_vectors. If not given (default), record all \(k\) left singular vectors.

svdsolverstr or callable

Strategy for computing the thin SVD of the states.

Options:

  • "dense" (default): Use scipy.linalg.svd() to compute the SVD. May be inefficient for very large state matrices.

  • "randomized": Compute an approximate SVD with a randomized approach via sklearn.utils.extmath.randomized_svd(). May be more efficient but less accurate for very large state matrices. NOTE: it is highly recommended to set max_vectors to limit the number of computed singular vectors. In this case, only max_vectors singular values are computed as well, meaning the cumulative and residual energies cannot be computed exactly.

  • callable: If this argument is a callable function, use it for the SVD computation. The signature must match scipy.linalg.svd(), i.e., U, s, Vh = svdsolver(states, **svdsolver_options)

weights(n, n) ndarray or (n,) ndarray None

Weight matrix \(\W\) or its diagonals. When provided, a weighted singular value decomposition of the states is used to ensure that the left singular vectors are orthogonal with respect to the weight matrix, i.e., \(\bfPhi\trp\W\bfPhi = \I\). If None (default), set \(\W\) to the identity.

namestr

Label for the state variable that this basis approximates.

svdsolver_optionsdict

Options to pass to the SVD solver.

Properties

cumulative_energy

Amount of singular value energy captured by the basis, \(\sum_{i=1}^r\sigma_i^2\big/\sum_{j=1}^k\sigma_j^2\).

entries

Entries of the basis matrix \(\Vr\in\RR^{n \times r}\).

full_state_dimension

Dimension \(n\) of the full state.

leftvecs

Leading left singular vectors of the training data.

max_vectors

Number of POD basis vectors stored in memory.

name

Label for the state variable that this basis approximates.

reduced_state_dimension

Dimension \(r\) of the reduced (compressed) state.

residual_energy

Amount of singular value energy not captured by the basis, \(\sum_{i=r+1}^k\sigma_i^2\big/\sum_{j=1}^k\sigma_j^2\).

rightvecs

Leading right singular vectors of the training data.

shape

Dimensions \((n, r)\) of the basis.

svdsolver

Strategy for computing the thin SVD of the states, either 'dense', 'randomized', or 'custom'.

svdsolver_options

Options to pass to the SVD solver.

svdvals

Normalized singular values of the training data.

weights

Weight matrix \(\W \in \RR^{n \times n}\).

Methods

compress

Map high-dimensional states to low-dimensional latent coordinates.

decompress

Map low-dimensional latent coordinates to high-dimensional states.

fit

Compute the POD basis entries by taking the SVD of the states.

from_svd

Initialize a PODBasis from the singular value decomposition of an \(n\times k\) matrix.

load

Load a basis from an HDF5 file.

plot1D

Plot the basis vectors over a one-dimensional domain.

plot_cumulative_energy

Plot the cumulative singular value energy.

plot_energy

Plot the normalized singular values and the cumulative and residual energies.

plot_residual_energy

Plot the residual singular value energy.

plot_svdval_decay

Plot the normalized singular value decay.

project

Project a high-dimensional state vector to the subset of the high-dimensional space that can be represented by the basis.

projection_error

Compute the error of the basis representation of a state or states.

save

Save the basis to an HDF5 file.

set_dimension

Set the reduced state dimension \(r\).

verify

Verify that compress() and decompress() are consistent in the sense that the range of decompress() is in the domain of compress() and that project() defines a projection operator, i.e., project(project(Q)) = project(Q).