flowtorch.analysis¶
flowtorch.analysis.svd¶
Classes and functions wrapping around torch.linalg.svd.
- class flowtorch.analysis.svd.SVD(data_matrix: Tensor, rank: Optional[int] = None, robust: Union[bool, dict] = False)[source]¶
Bases:
object
Compute and analyze the SVD of a data matrix.
- Parameters
U (pt.Tensor) – left singular vectors
s (pt.Tensor) – singular values
s_rel (pt.Tensor) – singular values normalized with their sum in percent
s_cum (pt.Tensor) – cumulative normalized singular values in percent
V (pt.Tensor) – right singular values
L (pt.Tensor) – low rank contribution to data matrix
S (pt.Tensor) – sparse contribution to data matrix
robust (Union[bool,dict]) – data_matrix is split into low rank and sparse contributions if True or if dictionary with options for Inexact ALM algorithm; the SVD is computed only on the low rank matrix; defaults to False
rank (int) – rank used for truncation
opt_rank (int) – optimal rank according to SVHT
required_memory (int) – memory required to store the truncation U, s, and V
Examples
>>> from flowtorch import DATASETS >>> from flowtorch.data import FOAMDataloader >>> from flowtorch.analysis import SVD >>> loader = FOAMDataloader(DATASETS["of_cavity_ascii"]) >>> data = loader.load_snapshot("p", loader.write_times[1:]) >>> svd = SVD(data, rank=100) >>> print(svd) SVD of a 400x5 data matrix Selected/optimal rank: 5/2 data type: torch.float32 (4b) truncated SVD size: 7.9297Kb >>> svd.s_rel tensor([9.9969e+01, 3.0860e-02, 3.0581e-04, 7.8097e-05, 3.2241e-05]) >>> svd.s_cum tensor([ 99.9687, 99.9996, 99.9999, 100.0000, 100.0000]) >>> svd.U.shape torch.Size([400, 5]) >>> svd = SVD(data, rank=100, robust=True) >>> svd.L.shape torch.Size([400, 5]) >>> svd = SVD(data, rank=100, robust={"sparsity" : 1.0, "verbose" : True, "max_iter" : 100}) >>> svd.S.shape torch.Size([400, 5])
- property L: Tensor¶
- property S: Tensor¶
- property U: Tensor¶
- property V: Tensor¶
- property opt_rank: int¶
- property rank: int¶
- reconstruct(rank: Optional[int] = None) Tensor [source]¶
Reconstruct the data matrix for a given rank.
- Parameters
rank (int, optional) – rank used to compute a truncated reconstruction
- Returns
reconstruction of the input data matrix
- Return type
pt.Tensor
- property required_memory: int¶
Compute the memory size in bytes of the truncated SVD.
- Returns
cumulative size of truncated U, s, and V tensors in bytes
- Return type
int
- property robust: Union[bool, dict]¶
- property s: Tensor¶
- property s_cum: Tensor¶
- property s_rel: Tensor¶
- flowtorch.analysis.svd.inexact_alm_matrix_complection(data_matrix: Tensor, sparsity: float = 1.0, tol: float = 1e-06, max_iter: int = 100, verbose: bool = False) Tuple[Tensor, Tensor] [source]¶
Split a data matrix in low rank and sparse contributions.
This function implements the inexact augmented Lagrange multiplier matrix completion algorithm to solve the principal component pursuit problem. The implementation is based on the Matlab code of Isabel Scherl (link), which is in turn based on the LRSLibrary.
- Parameters
data_matrix (pt.Tensor) – input data matrix; snapshots must be organized as column vectors
sparsity (float, optional) – factor to compute Lagrangian multiplyer for sparsity (typically named lambda); lower values lead to more agressive filtering
tol (float, optional) – tolerance for the normalized Frobenius norm of the difference between original data and the sum of low rank and sparse contributions; defaults to 1.0e-6
max_iter (int, optional) – maximum number of iteration before to give up, defaults to 100
verbose (bool, optional) – residual is printed for every iteration if True; defaults to False
- Returns
tuple holding the low rank and the sparse matrices
- Return type
Tuple[pt.Tensor, pt.Tensor]
flowtorch.analysis.dmd¶
Classes and functions to compute the dynamic mode decomposition (DMD) of a data matrix.
- class flowtorch.analysis.dmd.DMD(data_matrix: Tensor, dt: float, rank: Optional[int] = None, robust: Union[bool, dict] = False, unitary: bool = False, optimal: bool = False, tlsq=False)[source]¶
Bases:
object
Class computing the exact DMD of a data matrix.
Currently, no advanced mode selection algorithms are implemented. The mode amplitudes are computed using the first snapshot.
Examples
>>> from flowtorch import DATASETS >>> from flowtorch.data import FOAMDataloader >>> from flowtorch.analysis import DMD >>> path = DATASETS["of_cavity_binary"] >>> loader = FOAMDataloader(path) >>> data_matrix = loader.load_snapshot("p", loader.write_times) >>> dmd = DMD(data_matrix, dt=0.1, rank=3) >>> dmd.frequency tensor([0., 5., 0.]) >>> dmd.growth_rate tensor([-2.3842e-06, -4.2345e+01, -1.8552e+01]) >>> dmd.amplitude tensor([10.5635+0.j, -0.0616+0.j, -0.0537+0.j]) >>> dmd = DMD(data_matrix, dt=0.1, rank=3, robust=True) >>> dmd = DMD(data_matrix, dt=0.1, rank=3, robust={"tol": 1.0e-5, "verbose" : True})
- property amplitude: Tensor¶
- property dynamics: Tensor¶
- property eigvals: Tensor¶
- property eigvecs: Tensor¶
- property frequency: Tensor¶
- property growth_rate: Tensor¶
- property integral_contribution: Tensor¶
Integral contribution of individual modes according to J. Kou et al. 2017.
- property modes: Tensor¶
- property operator: Tensor¶
- partial_reconstruction(mode_indices: Set[int]) Tensor [source]¶
Reconstruct data matrix with limited number of modes.
- Parameters
mode_indices (Set[int]) – mode indices to keep
- Returns
reconstructed data matrix
- Return type
pt.Tensor
- property projection_error: Tensor¶
Compute the difference between Y and AX.
- Returns
projection error
- Return type
pt.Tensor
- property reconstruction: Tensor¶
Reconstruct an approximation of the training data.
- Returns
reconstructed training data
- Return type
pt.Tensor
- property reconstruction_error: Tensor¶
Compute the reconstruction error.
- Returns
difference between reconstruction and data matrix
- Return type
pt.Tensor
- property required_memory: int¶
Compute the memory size in bytes of the DMD.
- Returns
cumulative size of SVD, eigen values/vectors, and DMD modes in bytes
- Return type
int
- property tlsq_error: Tuple[Tensor, Tensor]¶
Compute the noise in X and Y.
- Returns
noise in X and Y
- Return type
Tuple[pt.Tensor, pt.Tensor]
- top_modes(n: int = 10, integral: bool = False) Tensor [source]¶
Get the indices of the first n most important modes.
Note that the conjugate complex modes for real data matrices are not filtered out.
- Parameters
n (int) – number of indices to return; defaults to 10
integral (bool, optional) – if True, the modes are sorted according to their integral contribution; defaults to False
- Returns
indices of top n modes sorted by amplitude or integral contribution
- Return type
pt.Tensor
flowtorch.analysis.psp_explorer¶
Module with classes and functions to explore and analyse iPSP data.
- class flowtorch.analysis.psp_explorer.PSPExplorer(file_path: str)[source]¶
Bases:
object
Explore iPSP data interactively.
- interact(zone: str, field_name: str, times: list, mask: bool = True, width: int = 1024, every: int = 5, cmin: float = -2, cmax: float = 0.5) Figure [source]¶
Show selected snapshots of an iPSP file as interactive surface plot.
- Parameters
zone (str) – name of the zone to display
field_name (str) – name of the field to display
times (list) – list of times instances to show
mask (bool, optional) – if True, the data is masked using the mask provided in the iPSP dataset, defaults to True
width (int, optional) – width of the plot, defaults to 1024
every (int, optional) – plot only every nth point, defaults to 5
cmin (float, optional) – lower bound for color range, defaults to -2
cmax (float, optional) – upper bound for color range, defaults to 0.5
- Returns
interactive iPSP surface plot
- Return type
go.Figure
- property loader: PSPDataloader¶
Get the PSPdataloader instance used to access data.
This propertie allows accessing the loader’s metadata via the explorer without having to create a new loader instance.
- Returns
PSP dataloader instance
- Return type
- mean(zone: str, field_name: str, times: list, mask: bool = True, width: int = 1024, every=5, cmin: float = -2, cmax: float = 0.5) Figure [source]¶
Show temporal mean of an iPSP file as interactive surface plot.
- Parameters
zone (str) – name of the zone to display
field_name (str) – name of the field to display
times (list) – snapshot times over which to compute the mean
mask (bool, optional) – if True, the data is masked using the mask provided in the iPSP dataset, defaults to True
width (int, optional) – width of the plot, defaults to 1024
every (int, optional) – use only every nth point for plotting, defaults to 5
cmin (float, optional) – lower bound for color range, defaults to -2
cmax (float, optional) – upper bound for color range, defaults to 0.5
- Returns
interactive iPSP surface plot showing the temporal mean
- Return type
go.Figure
- std(zone: str, field_name: str, times: list, mask: bool = True, width: int = 1024, every=5, cmin: float = 0, cmax: float = 0.5) Figure [source]¶
Show temporal standard deviation as interactive surface plot.
- Parameters
zone (str) – name of the zone to display
field_name (str) – name of the field to display
times (list) – snapshot times over which to compute the standard deviation
mask (bool, optional) – if True, the data is masked using the mask provided in the iPSP dataset, defaults to True
width (int, optional) – width of the plot, defaults to 1024
every (int, optional) – use only every nth point for plotting, defaults to 5
cmin (float, optional) – lower bound for color range, defaults to -2
cmax (float, optional) – upper bound for color range, defaults to 0.5
- Returns
interactive iPSP surface plot
- Return type
go.Figure