romtools.vector_space
This module defines the API to work with a vector subspace. A vector subspace is foundational to reduced-order models. In a ROM, a high-dimensional state is restricted to live within a low-dimensional vector space, known as a trial space. Mathematically, given a "FOM" vector $\mathbf{u} \in \mathbb{R}^{N_{\mathrm{vars}} N_{\mathrm{x}}}$, we can write $$\mathbf{u} \approx \tilde{\mathbf{u}} \in \mathcal{V} + \mathbf{u}_{\mathrm{shift}}$$ where
- $\mathcal{V}$ with $\text{dim}(\mathcal{V}) = K \le N_{\mathrm{vars}} N_{\mathrm{x}}$ is the trial space
- $N_{\mathrm{vars}}$ is the number of PDE variables (e.g., 5 for the compressible Navier-Stokes equations in 3D)
- $N_{\mathrm{x}}$ is the number of spatial DOFs
Formally, we can describe this low-dimensional representation with a basis and an affine offset, $$\tilde{\mathbf{u}} = \boldsymbol \Phi \hat{\mathbf{u}} + \mathbf{u}_{\mathrm{shift}}$$ where $\boldsymbol \Phi \in \mathbb{R}^{ N_{\mathrm{vars}} N_{\mathrm{x}} \times K}$ is the basis matrix ($K$ is the number of basis), $\hat{\mathbf{u}} \in \mathbb{R}^{K}$ are the reduced, or generalized coordinates, $\mathbf{u}_{\mathrm{shift}} \in \mathbb{R}^{ N_{\mathrm{vars}} N_{\mathrm{x}}}$ is the shift vector (or affine offset), and, by definition, $\mathcal{V} \equiv \mathrm{range}(\boldsymbol \Phi)$.
The VectorSpace
abstract class defined below encapsulates the information of an affine vector space, $\mathcal{V}$,
by virtue of providing access to a basis matrix, a shift vector, and the dimensionality of the vector space,
while decoupling this representation from how it is computed.
We rely on a tensor representation!
Our representation of the basis and the affine offset for a vector space is based on tensors $$\mathcal{\Phi} \in \mathbb{R}^{ N_{\mathrm{vars}} \times N_{\mathrm{x}} \times K},$$ $$\mathcal{u}_{\mathrm{shift}} \in \mathbb{R}^{ N_{\mathrm{vars}} \times N_{\mathrm{x}}}.$$ Internally, we remark that all tensors are reshaped into 2D matrices, e.g., when performing SVD.
Content
We currently provide the following concrete classes:
DictionaryVectorSpace
: construct a vector space from a matrix without truncation.VectorSpaceFromPOD
: construct a vector subspace computed via SVD.
which derive from the abstract class VectorSpace
.
API
1# 2# ************************************************************************ 3# 4# ROM Tools and Workflows 5# Copyright 2019 National Technology & Engineering Solutions of Sandia,LLC 6# (NTESS) 7# 8# Under the terms of Contract DE-NA0003525 with NTESS, the 9# U.S. Government retains certain rights in this software. 10# 11# ROM Tools and Workflows is licensed under BSD-3-Clause terms of use: 12# 13# Redistribution and use in source and binary forms, with or without 14# modification, are permitted provided that the following conditions 15# are met: 16# 17# 1. Redistributions of source code must retain the above copyright 18# notice, this list of conditions and the following disclaimer. 19# 20# 2. Redistributions in binary form must reproduce the above copyright 21# notice, this list of conditions and the following disclaimer in the 22# documentation and/or other materials provided with the distribution. 23# 24# 3. Neither the name of the copyright holder nor the names of its 25# contributors may be used to endorse or promote products derived 26# from this software without specific prior written permission. 27# 28# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS 29# "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT 30# LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS 31# FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE 32# COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, 33# INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES 34# (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR 35# SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) 36# HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, 37# STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING 38# IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE 39# POSSIBILITY OF SUCH DAMAGE. 40# 41# Questions? Contact Eric Parish (ejparis@sandia.gov) 42# 43# ************************************************************************ 44# 45 46''' 47This module defines the API to work with a vector subspace. 48A vector subspace is foundational to reduced-order models. 49In a ROM, a high-dimensional state is restricted to live within a low-dimensional vector space, known as a trial space. 50Mathematically, given a "FOM" vector $\\mathbf{u} \\in \\mathbb{R}^{N_{\\mathrm{vars}} N_{\\mathrm{x}}}$, we can write 51$$\\mathbf{u} \\approx \\tilde{\\mathbf{u}} \\in \\mathcal{V} + \\mathbf{u}_{\\mathrm{shift}}$$ 52where 53- $\\mathcal{V}$ with $\\text{dim}(\\mathcal{V}) = K \\le N_{\\mathrm{vars}} N_{\\mathrm{x}}$ is the trial space 54- $N_{\\mathrm{vars}}$ is the number of PDE variables (e.g., 5 for the compressible Navier-Stokes equations in 3D) 55- $N_{\\mathrm{x}}$ is the number of spatial DOFs 56 57Formally, we can describe this low-dimensional representation with a basis and an affine offset, 58$$\\tilde{\\mathbf{u}} = \\boldsymbol \\Phi \\hat{\\mathbf{u}} + \\mathbf{u}_{\\mathrm{shift}}$$ 59where $\\boldsymbol \\Phi \\in \\mathbb{R}^{ N_{\\mathrm{vars}} N_{\\mathrm{x}} \\times K}$ is the basis matrix 60($K$ is the number of basis), $\\hat{\\mathbf{u}} \\in \\mathbb{R}^{K}$ are the reduced, or generalized coordinates, 61$\\mathbf{u}_{\\mathrm{shift}} \\in \\mathbb{R}^{ N_{\\mathrm{vars}} N_{\\mathrm{x}}}$ is the shift vector (or affine offset), 62and, by definition, $\\mathcal{V} \\equiv \\mathrm{range}(\\boldsymbol \\Phi)$. 63 64The `VectorSpace` abstract class defined below encapsulates the information of an affine vector space, $\\mathcal{V}$, 65by virtue of providing access to a basis matrix, a shift vector, and the dimensionality of the vector space, 66while decoupling this representation from *how* it is computed. 67 68####We rely on a tensor representation! 69 70Our representation of the basis and the affine offset for a vector space is based on tensors 71$$\\mathcal{\\Phi} \\in \\mathbb{R}^{ N_{\\mathrm{vars}} \\times N_{\\mathrm{x}} \\times K},$$ 72$$\\mathcal{u}_{\\mathrm{shift}} \\in \\mathbb{R}^{ N_{\\mathrm{vars}} \\times N_{\\mathrm{x}}}.$$ 73Internally, we remark that all tensors are reshaped into 2D matrices, e.g., when performing SVD. 74 75###**Content** 76 77We currently provide the following concrete classes: 78 79- `DictionaryVectorSpace`: construct a vector space from a matrix without truncation. 80 81- `VectorSpaceFromPOD`: construct a vector subspace computed via SVD. 82 83which derive from the abstract class `VectorSpace`. 84 85--- 86##**API** 87''' 88 89from typing import Tuple, Protocol, Callable 90import numpy as np 91from romtools.vector_space.utils.truncater import LeftSingularVectorTruncater, NoOpTruncater 92from romtools.vector_space.utils.shifter import Shifter, create_noop_shifter 93from romtools.vector_space.utils.scaler import Scaler, NoOpScaler 94from romtools.vector_space.utils.orthogonalizer import Orthogonalizer, NoOpOrthogonalizer 95 96 97class VectorSpace(Protocol): 98 ''' 99 Abstract base class for vector space implementations. 100 101 Methods: 102 ''' 103 104 def extents(self) -> Tuple[int, int, int]: 105 ''' 106 Retrieves the dimension of the vector space 107 108 Returns: 109 A Tuple with the the dimensions of the vector space (n_var,nx,K). 110 ''' 111 ... 112 113 def get_shift_vector(self) -> np.ndarray: 114 ''' 115 Retrieves the shift vector of the vector space. 116 117 Returns: 118 `np.ndarray`: The shift vector in tensorm form. 119 ''' 120 ... 121 122 def get_basis(self) -> np.ndarray: 123 ''' 124 Retrieves the basis vectors of the vector space. 125 126 Returns: 127 `np.ndarray`: The basis of the vector space in tensor form. 128 ''' 129 ... 130 131 132class DictionaryVectorSpace(): 133 ''' 134 Reduced basis vector space (no truncation). 135 136 This class conforms to `VectorSpace` protocol. 137 138 Given a snapshot matrix $\\mathbf{S}$, we set the basis to be 139 140 $$\\boldsymbol \\Phi = \\mathrm{orthogonalize}(\\mathbf{S} - \\mathbf{u}_{\\mathrm{shift}})$$ 141 142 where the orthogonalization and shifts are defined by their 143 respective classes 144 ''' 145 146 def __init__(self, 147 snapshots, 148 shifter: Shifter = None, 149 orthogonalizer: Orthogonalizer = NoOpOrthogonalizer()) -> None: 150 ''' 151 Constructor. 152 153 Args: 154 snapshots (np.ndarray): Snapshot data in tensor form 155 $\\in \\mathbb{R}^{ N_{\\mathrm{vars}} \\times N_{\\mathrm{x}} \\times N_{samples}}$ 156 shifter: Class that shifts the basis. 157 orthogonalizer: Class that orthogonalizes the basis. 158 159 This constructor initializes a vector space by performing basis 160 manipulation operations on the provided snapshot data. 161 ''' 162 # Create noop shifter if not provided 163 if shifter is None: 164 shifter = create_noop_shifter(snapshots) 165 166 # compute basis 167 n_var = snapshots.shape[0] 168 shifter.apply_shift(snapshots) 169 snapshot_matrix = _tensor_to_matrix(snapshots) 170 self.__basis = orthogonalizer.orthogonalize(snapshot_matrix) 171 self.__basis = _matrix_to_tensor(n_var, self.__basis) 172 self.__shift_vector = shifter.get_shift_vector() 173 174 def get_shift_vector(self) -> np.ndarray: 175 ''' 176 Concrete implementation of `VectorSpace.get_shift_vector()` 177 ''' 178 return self.__shift_vector 179 180 def get_basis(self) -> np.ndarray: 181 ''' 182 Concrete implementation of `VectorSpace.get_basis()` 183 ''' 184 return self.__basis 185 186 def extents(self) -> Tuple[int, int, int]: 187 return self.__basis.shape 188 189 190class VectorSpaceFromPOD(): 191 ''' 192 POD vector space (constructed via SVD). 193 194 This class conforms to `VectorSpace` protocol. 195 196 Given a snapshot matrix $\\mathbf{S}$, we compute the basis $\\boldsymbol \\Phi$ as 197 198 199 $$\\boldsymbol U = \\mathrm{SVD}(\\mathrm{prescale}(\\mathbf{S} - \\mathbf{u}_{\\mathrm{shift}}))$$ 200 $$\\boldsymbol \\Phi = \\mathrm{orthogonalize}(\\mathrm{postscale}(\\mathrm{truncate}( \\boldsymbol U )))$$ 201 202 where $\\boldsymbol U$ are the left singular vectors and the 203 orthogonalization, truncation, scaling, and shifts are defined by their respective classes. 204 205 For truncation, we enable truncation based on a fixed dimension or the 206 decay of singular values; please refer to the documentation for the truncater. 207 ''' 208 209 def __init__(self, 210 snapshots, 211 truncater: LeftSingularVectorTruncater = NoOpTruncater(), 212 shifter: Shifter = None, 213 orthogonalizer: Orthogonalizer = NoOpOrthogonalizer(), 214 scaler: Scaler = NoOpScaler(), 215 svdFnc: Callable = None) -> None: 216 ''' 217 Constructor. 218 219 Args: 220 snapshots (np.ndarray): Snapshot data in tensor form 221 $\\in \\mathbb{R}^{ N_{\\mathrm{vars}} \\times N_{\\mathrm{x}} \\times N_{samples}}$ 222 truncater (Truncater): Concrete implementation for truncating the basis. 223 shifter (Shifter): Concrete implementation responsible for shifting the basis. 224 orthogonalizer (Orthogonalizer): Concrete implementation that orthogonalizes the basis. 225 scaler: Concrete implementation that scales the basis. 226 svdFnc: a callable to use for computing the SVD on the snapshots data. 227 IMPORTANT: must conform to the API of 228 [np.linalg.svd](https://numpy.org/doc/stable/reference/generated/numpy.linalg.svd.html#numpy-linalg-svd). 229 If `None`, internally we use `np.linalg.svd`. 230 Note: this is useful when you want to use a custom svd, for 231 example when your snapshots are distributed with MPI, or 232 maybe you have a fancy svd function that you can use. 233 234 This constructor initializes a POD vector space by performing SVD on the provided snapshot 235 data and applying various basis manipulation operations, including truncation, shifting, scaling, 236 and orthogonalization. 237 ''' 238 if shifter is None: 239 shifter = create_noop_shifter(snapshots) 240 n_var = snapshots.shape[0] 241 shifter.apply_shift(snapshots) 242 scaler.pre_scale(snapshots) 243 snapshot_matrix = _tensor_to_matrix(snapshots) 244 svd_picked = np.linalg.svd if svdFnc is None else svdFnc 245 lsv, svals, _ = svd_picked(snapshot_matrix, full_matrices=False, 246 compute_uv=True, hermitian=False) 247 248 self.__basis = truncater.truncate(lsv, svals) 249 self.__basis = _matrix_to_tensor(n_var, self.__basis) 250 scaler.post_scale(self.__basis) 251 self.__basis = _tensor_to_matrix(self.__basis) 252 self.__basis = orthogonalizer.orthogonalize(self.__basis) 253 self.__basis = _matrix_to_tensor(n_var, self.__basis) 254 self.__shift_vector = shifter.get_shift_vector() 255 256 def get_shift_vector(self) -> np.ndarray: 257 ''' 258 Concrete implementation of `VectorSpace.get_shift_vector()` 259 ''' 260 return self.__shift_vector 261 262 def get_basis(self) -> np.ndarray: 263 ''' 264 Concrete implementation of `VectorSpace.get_basis()` 265 ''' 266 return self.__basis 267 268 def extents(self) -> Tuple[int, int, int]: 269 return self.__basis.shape 270 271 272def _tensor_to_matrix(tensor_input: np.ndarray) -> np.ndarray: 273 ''' 274 Converts a tensor with shape $[N, M, P]$ to a matrix representation 275 in which the first two dimension are collapsed $[N M, P]$. 276 ''' 277 output_tensor = tensor_input.reshape(tensor_input.shape[0]*tensor_input.shape[1], 278 tensor_input.shape[2]) 279 return output_tensor 280 281 282def _matrix_to_tensor(n_var: int, matrix_input: np.ndarray) -> np.ndarray: 283 ''' 284 Inverse operation of `_tensor_to_matrix` 285 ''' 286 d1 = int(matrix_input.shape[0] / n_var) 287 d2 = matrix_input.shape[1] 288 output_matrix = matrix_input.reshape(n_var, d1, d2) 289 return output_matrix
98class VectorSpace(Protocol): 99 ''' 100 Abstract base class for vector space implementations. 101 102 Methods: 103 ''' 104 105 def extents(self) -> Tuple[int, int, int]: 106 ''' 107 Retrieves the dimension of the vector space 108 109 Returns: 110 A Tuple with the the dimensions of the vector space (n_var,nx,K). 111 ''' 112 ... 113 114 def get_shift_vector(self) -> np.ndarray: 115 ''' 116 Retrieves the shift vector of the vector space. 117 118 Returns: 119 `np.ndarray`: The shift vector in tensorm form. 120 ''' 121 ... 122 123 def get_basis(self) -> np.ndarray: 124 ''' 125 Retrieves the basis vectors of the vector space. 126 127 Returns: 128 `np.ndarray`: The basis of the vector space in tensor form. 129 ''' 130 ...
Abstract base class for vector space implementations.
Methods:
1771def _no_init_or_replace_init(self, *args, **kwargs): 1772 cls = type(self) 1773 1774 if cls._is_protocol: 1775 raise TypeError('Protocols cannot be instantiated') 1776 1777 # Already using a custom `__init__`. No need to calculate correct 1778 # `__init__` to call. This can lead to RecursionError. See bpo-45121. 1779 if cls.__init__ is not _no_init_or_replace_init: 1780 return 1781 1782 # Initially, `__init__` of a protocol subclass is set to `_no_init_or_replace_init`. 1783 # The first instantiation of the subclass will call `_no_init_or_replace_init` which 1784 # searches for a proper new `__init__` in the MRO. The new `__init__` 1785 # replaces the subclass' old `__init__` (ie `_no_init_or_replace_init`). Subsequent 1786 # instantiation of the protocol subclass will thus use the new 1787 # `__init__` and no longer call `_no_init_or_replace_init`. 1788 for base in cls.__mro__: 1789 init = base.__dict__.get('__init__', _no_init_or_replace_init) 1790 if init is not _no_init_or_replace_init: 1791 cls.__init__ = init 1792 break 1793 else: 1794 # should not happen 1795 cls.__init__ = object.__init__ 1796 1797 cls.__init__(self, *args, **kwargs)
105 def extents(self) -> Tuple[int, int, int]: 106 ''' 107 Retrieves the dimension of the vector space 108 109 Returns: 110 A Tuple with the the dimensions of the vector space (n_var,nx,K). 111 ''' 112 ...
Retrieves the dimension of the vector space
Returns:
A Tuple with the the dimensions of the vector space (n_var,nx,K).
114 def get_shift_vector(self) -> np.ndarray: 115 ''' 116 Retrieves the shift vector of the vector space. 117 118 Returns: 119 `np.ndarray`: The shift vector in tensorm form. 120 ''' 121 ...
Retrieves the shift vector of the vector space.
Returns:
np.ndarray
: The shift vector in tensorm form.
123 def get_basis(self) -> np.ndarray: 124 ''' 125 Retrieves the basis vectors of the vector space. 126 127 Returns: 128 `np.ndarray`: The basis of the vector space in tensor form. 129 ''' 130 ...
Retrieves the basis vectors of the vector space.
Returns:
np.ndarray
: The basis of the vector space in tensor form.
133class DictionaryVectorSpace(): 134 ''' 135 Reduced basis vector space (no truncation). 136 137 This class conforms to `VectorSpace` protocol. 138 139 Given a snapshot matrix $\\mathbf{S}$, we set the basis to be 140 141 $$\\boldsymbol \\Phi = \\mathrm{orthogonalize}(\\mathbf{S} - \\mathbf{u}_{\\mathrm{shift}})$$ 142 143 where the orthogonalization and shifts are defined by their 144 respective classes 145 ''' 146 147 def __init__(self, 148 snapshots, 149 shifter: Shifter = None, 150 orthogonalizer: Orthogonalizer = NoOpOrthogonalizer()) -> None: 151 ''' 152 Constructor. 153 154 Args: 155 snapshots (np.ndarray): Snapshot data in tensor form 156 $\\in \\mathbb{R}^{ N_{\\mathrm{vars}} \\times N_{\\mathrm{x}} \\times N_{samples}}$ 157 shifter: Class that shifts the basis. 158 orthogonalizer: Class that orthogonalizes the basis. 159 160 This constructor initializes a vector space by performing basis 161 manipulation operations on the provided snapshot data. 162 ''' 163 # Create noop shifter if not provided 164 if shifter is None: 165 shifter = create_noop_shifter(snapshots) 166 167 # compute basis 168 n_var = snapshots.shape[0] 169 shifter.apply_shift(snapshots) 170 snapshot_matrix = _tensor_to_matrix(snapshots) 171 self.__basis = orthogonalizer.orthogonalize(snapshot_matrix) 172 self.__basis = _matrix_to_tensor(n_var, self.__basis) 173 self.__shift_vector = shifter.get_shift_vector() 174 175 def get_shift_vector(self) -> np.ndarray: 176 ''' 177 Concrete implementation of `VectorSpace.get_shift_vector()` 178 ''' 179 return self.__shift_vector 180 181 def get_basis(self) -> np.ndarray: 182 ''' 183 Concrete implementation of `VectorSpace.get_basis()` 184 ''' 185 return self.__basis 186 187 def extents(self) -> Tuple[int, int, int]: 188 return self.__basis.shape
Reduced basis vector space (no truncation).
This class conforms to VectorSpace
protocol.
Given a snapshot matrix $\mathbf{S}$, we set the basis to be
$$\boldsymbol \Phi = \mathrm{orthogonalize}(\mathbf{S} - \mathbf{u}_{\mathrm{shift}})$$
where the orthogonalization and shifts are defined by their respective classes
147 def __init__(self, 148 snapshots, 149 shifter: Shifter = None, 150 orthogonalizer: Orthogonalizer = NoOpOrthogonalizer()) -> None: 151 ''' 152 Constructor. 153 154 Args: 155 snapshots (np.ndarray): Snapshot data in tensor form 156 $\\in \\mathbb{R}^{ N_{\\mathrm{vars}} \\times N_{\\mathrm{x}} \\times N_{samples}}$ 157 shifter: Class that shifts the basis. 158 orthogonalizer: Class that orthogonalizes the basis. 159 160 This constructor initializes a vector space by performing basis 161 manipulation operations on the provided snapshot data. 162 ''' 163 # Create noop shifter if not provided 164 if shifter is None: 165 shifter = create_noop_shifter(snapshots) 166 167 # compute basis 168 n_var = snapshots.shape[0] 169 shifter.apply_shift(snapshots) 170 snapshot_matrix = _tensor_to_matrix(snapshots) 171 self.__basis = orthogonalizer.orthogonalize(snapshot_matrix) 172 self.__basis = _matrix_to_tensor(n_var, self.__basis) 173 self.__shift_vector = shifter.get_shift_vector()
Constructor.
Arguments:
- snapshots (np.ndarray): Snapshot data in tensor form $\in \mathbb{R}^{ N_{\mathrm{vars}} \times N_{\mathrm{x}} \times N_{samples}}$
- shifter: Class that shifts the basis.
- orthogonalizer: Class that orthogonalizes the basis.
This constructor initializes a vector space by performing basis manipulation operations on the provided snapshot data.
175 def get_shift_vector(self) -> np.ndarray: 176 ''' 177 Concrete implementation of `VectorSpace.get_shift_vector()` 178 ''' 179 return self.__shift_vector
Concrete implementation of VectorSpace.get_shift_vector()
191class VectorSpaceFromPOD(): 192 ''' 193 POD vector space (constructed via SVD). 194 195 This class conforms to `VectorSpace` protocol. 196 197 Given a snapshot matrix $\\mathbf{S}$, we compute the basis $\\boldsymbol \\Phi$ as 198 199 200 $$\\boldsymbol U = \\mathrm{SVD}(\\mathrm{prescale}(\\mathbf{S} - \\mathbf{u}_{\\mathrm{shift}}))$$ 201 $$\\boldsymbol \\Phi = \\mathrm{orthogonalize}(\\mathrm{postscale}(\\mathrm{truncate}( \\boldsymbol U )))$$ 202 203 where $\\boldsymbol U$ are the left singular vectors and the 204 orthogonalization, truncation, scaling, and shifts are defined by their respective classes. 205 206 For truncation, we enable truncation based on a fixed dimension or the 207 decay of singular values; please refer to the documentation for the truncater. 208 ''' 209 210 def __init__(self, 211 snapshots, 212 truncater: LeftSingularVectorTruncater = NoOpTruncater(), 213 shifter: Shifter = None, 214 orthogonalizer: Orthogonalizer = NoOpOrthogonalizer(), 215 scaler: Scaler = NoOpScaler(), 216 svdFnc: Callable = None) -> None: 217 ''' 218 Constructor. 219 220 Args: 221 snapshots (np.ndarray): Snapshot data in tensor form 222 $\\in \\mathbb{R}^{ N_{\\mathrm{vars}} \\times N_{\\mathrm{x}} \\times N_{samples}}$ 223 truncater (Truncater): Concrete implementation for truncating the basis. 224 shifter (Shifter): Concrete implementation responsible for shifting the basis. 225 orthogonalizer (Orthogonalizer): Concrete implementation that orthogonalizes the basis. 226 scaler: Concrete implementation that scales the basis. 227 svdFnc: a callable to use for computing the SVD on the snapshots data. 228 IMPORTANT: must conform to the API of 229 [np.linalg.svd](https://numpy.org/doc/stable/reference/generated/numpy.linalg.svd.html#numpy-linalg-svd). 230 If `None`, internally we use `np.linalg.svd`. 231 Note: this is useful when you want to use a custom svd, for 232 example when your snapshots are distributed with MPI, or 233 maybe you have a fancy svd function that you can use. 234 235 This constructor initializes a POD vector space by performing SVD on the provided snapshot 236 data and applying various basis manipulation operations, including truncation, shifting, scaling, 237 and orthogonalization. 238 ''' 239 if shifter is None: 240 shifter = create_noop_shifter(snapshots) 241 n_var = snapshots.shape[0] 242 shifter.apply_shift(snapshots) 243 scaler.pre_scale(snapshots) 244 snapshot_matrix = _tensor_to_matrix(snapshots) 245 svd_picked = np.linalg.svd if svdFnc is None else svdFnc 246 lsv, svals, _ = svd_picked(snapshot_matrix, full_matrices=False, 247 compute_uv=True, hermitian=False) 248 249 self.__basis = truncater.truncate(lsv, svals) 250 self.__basis = _matrix_to_tensor(n_var, self.__basis) 251 scaler.post_scale(self.__basis) 252 self.__basis = _tensor_to_matrix(self.__basis) 253 self.__basis = orthogonalizer.orthogonalize(self.__basis) 254 self.__basis = _matrix_to_tensor(n_var, self.__basis) 255 self.__shift_vector = shifter.get_shift_vector() 256 257 def get_shift_vector(self) -> np.ndarray: 258 ''' 259 Concrete implementation of `VectorSpace.get_shift_vector()` 260 ''' 261 return self.__shift_vector 262 263 def get_basis(self) -> np.ndarray: 264 ''' 265 Concrete implementation of `VectorSpace.get_basis()` 266 ''' 267 return self.__basis 268 269 def extents(self) -> Tuple[int, int, int]: 270 return self.__basis.shape
POD vector space (constructed via SVD).
This class conforms to VectorSpace
protocol.
Given a snapshot matrix $\mathbf{S}$, we compute the basis $\boldsymbol \Phi$ as
$$\boldsymbol U = \mathrm{SVD}(\mathrm{prescale}(\mathbf{S} - \mathbf{u}_{\mathrm{shift}}))$$ $$\boldsymbol \Phi = \mathrm{orthogonalize}(\mathrm{postscale}(\mathrm{truncate}( \boldsymbol U )))$$
where $\boldsymbol U$ are the left singular vectors and the orthogonalization, truncation, scaling, and shifts are defined by their respective classes.
For truncation, we enable truncation based on a fixed dimension or the decay of singular values; please refer to the documentation for the truncater.
210 def __init__(self, 211 snapshots, 212 truncater: LeftSingularVectorTruncater = NoOpTruncater(), 213 shifter: Shifter = None, 214 orthogonalizer: Orthogonalizer = NoOpOrthogonalizer(), 215 scaler: Scaler = NoOpScaler(), 216 svdFnc: Callable = None) -> None: 217 ''' 218 Constructor. 219 220 Args: 221 snapshots (np.ndarray): Snapshot data in tensor form 222 $\\in \\mathbb{R}^{ N_{\\mathrm{vars}} \\times N_{\\mathrm{x}} \\times N_{samples}}$ 223 truncater (Truncater): Concrete implementation for truncating the basis. 224 shifter (Shifter): Concrete implementation responsible for shifting the basis. 225 orthogonalizer (Orthogonalizer): Concrete implementation that orthogonalizes the basis. 226 scaler: Concrete implementation that scales the basis. 227 svdFnc: a callable to use for computing the SVD on the snapshots data. 228 IMPORTANT: must conform to the API of 229 [np.linalg.svd](https://numpy.org/doc/stable/reference/generated/numpy.linalg.svd.html#numpy-linalg-svd). 230 If `None`, internally we use `np.linalg.svd`. 231 Note: this is useful when you want to use a custom svd, for 232 example when your snapshots are distributed with MPI, or 233 maybe you have a fancy svd function that you can use. 234 235 This constructor initializes a POD vector space by performing SVD on the provided snapshot 236 data and applying various basis manipulation operations, including truncation, shifting, scaling, 237 and orthogonalization. 238 ''' 239 if shifter is None: 240 shifter = create_noop_shifter(snapshots) 241 n_var = snapshots.shape[0] 242 shifter.apply_shift(snapshots) 243 scaler.pre_scale(snapshots) 244 snapshot_matrix = _tensor_to_matrix(snapshots) 245 svd_picked = np.linalg.svd if svdFnc is None else svdFnc 246 lsv, svals, _ = svd_picked(snapshot_matrix, full_matrices=False, 247 compute_uv=True, hermitian=False) 248 249 self.__basis = truncater.truncate(lsv, svals) 250 self.__basis = _matrix_to_tensor(n_var, self.__basis) 251 scaler.post_scale(self.__basis) 252 self.__basis = _tensor_to_matrix(self.__basis) 253 self.__basis = orthogonalizer.orthogonalize(self.__basis) 254 self.__basis = _matrix_to_tensor(n_var, self.__basis) 255 self.__shift_vector = shifter.get_shift_vector()
Constructor.
Arguments:
- snapshots (np.ndarray): Snapshot data in tensor form $\in \mathbb{R}^{ N_{\mathrm{vars}} \times N_{\mathrm{x}} \times N_{samples}}$
- truncater (Truncater): Concrete implementation for truncating the basis.
- shifter (Shifter): Concrete implementation responsible for shifting the basis.
- orthogonalizer (Orthogonalizer): Concrete implementation that orthogonalizes the basis.
- scaler: Concrete implementation that scales the basis.
- svdFnc: a callable to use for computing the SVD on the snapshots data.
IMPORTANT: must conform to the API of
np.linalg.svd.
If
None
, internally we usenp.linalg.svd
. Note: this is useful when you want to use a custom svd, for example when your snapshots are distributed with MPI, or maybe you have a fancy svd function that you can use.
This constructor initializes a POD vector space by performing SVD on the provided snapshot data and applying various basis manipulation operations, including truncation, shifting, scaling, and orthogonalization.
257 def get_shift_vector(self) -> np.ndarray: 258 ''' 259 Concrete implementation of `VectorSpace.get_shift_vector()` 260 ''' 261 return self.__shift_vector
Concrete implementation of VectorSpace.get_shift_vector()