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 offit()
.The POD basis entries matrix \(\Vr = \bfPhi_{:,:r}\in\RR^{n\times r}\) always has orthonormal columns, 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 columns of the basis entries are also the dominant eigenvectors of \(\Q\trp\W\Q\) and can be computed through eigendecomposition by setting
svdsolver="eigh"
.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
, orprojection_error
. Once the basis entries are set by callingfit()
, the reduced state dimension \(r\) can be updated by callingset_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()
, thereduced_state_dimension
can be increased up tomax_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): Usescipy.linalg.svd()
to compute the SVD. May be inefficient for very large state matrices."randomized"
: Compute an approximate SVD with a randomized approach viasklearn.utils.extmath.randomized_svd()
. May be more efficient but less accurate for very large state matrices. NOTE: it is highly recommended to setmax_vectors
to limit the number of computed singular vectors. In this case, onlymax_vectors
singular values are computed as well, meaning the cumulative and residual energies cannot be computed exactly."method-of-snapshots"
or"eigh"
: Compute the basis through a symmetric eigenvalue decomposition, rather than through the SVD, viascipy.linalg.eigh()
. This is how POD was computed when it was orginally introduced. If the state dimension is larger than the number of snapshots, this method is much more efficient than the SVD. Moreover, non-Euclidean inner products (seeweights
) are handled much more efficiently this way than with an SVD-based approach. NOTE: in this case, an additional keyword argumentminthresh
defines a threshold at which small eigenvalues are truncated.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}\). Also accessible via indexing (
basis[:]
).
- 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. The
reduced_state_dimension
may be increased up tomax_vectors
.
- 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, if available.
- 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:Map high-dimensional states to low-dimensional latent coordinates.
Map low-dimensional latent coordinates to high-dimensional states.
Compute the POD basis entries by taking the SVD of the states.
Initialize a
PODBasis
from the singular value decomposition of an \(n\times k\) matrix.Load a basis from an HDF5 file.
Plot the basis vectors over a one-dimensional domain.
Plot the cumulative singular value energy and a function of the basis size.
Plot the normalized singular values and the cumulative and residual energies.
Plot the relative projection error of the training snapshots as a function of the basis size.
Plot the residual singular value energy as a function of the basis size.
Plot the normalized singular value decay.
Project a high-dimensional state vector to the subset of the high-dimensional space that can be represented by the basis.
Compute the error of the basis representation of a state or states.
Save the basis to an HDF5 file.
Set the reduced state dimension \(r\).
Verify that
compress()
anddecompress()
are consistent in the sense that the range ofdecompress()
is in the domain ofcompress()
and thatproject()
defines a projection operator, i.e.,project(project(Q)) = project(Q)
.