Towards customization: accessing MODULO internal modules¶
The methods for each decomposition (POD, DMD, mPOD, etc) are made of different sub-methods. It is possible to call each of these separately, hence allowing the user to customize the various steps if needed.
Functions for the POD¶
For example, the method compute_POD_K acts in three steps: 1. compute K, 2. diagonalize it to compute the Psi’s and 3 project to compute the Phi’s.
These are handled by three methods, with use the functions in modulo.core:
- modulo_vki.core._k_matrix.CorrelationMatrix(N_T, N_PARTITIONS=1, MEMORY_SAVING=False, FOLDER_OUT='./', SAVE_K=False, D=None, weights=array([], dtype=float64))[source]¶
This method computes the temporal correlation matrix, given a data matrix as input. It’s possible to use memory saving then splitting the computing in different tranches if computationally heavy. If D has been computed using MODULO then the dimension dim_col and N_PARTITIONS is automatically loaded
- Parameters:
N_T – int. Number of temporal snapshots
D – np.array. Data matrix
SAVE_K – bool. If SAVE_K=True, the matrix K is saved on disk. If the MEMORY_SAVING feature is active, this is done by default.
MEMORY_SAVING – bool. If MEMORY_SAVING = True, the computation of the correlation matrix is done by steps. It requires the data matrix to be partitioned, following algorithm in MODULO._data_processing.
FOLDER_OUT – str. Folder in which the temporal correlation matrix will be stored
N_PARTITIONS – int. Number of partitions to be read in computing the correlation matrix. If _data_processing is used to partition the data matrix, this is inherited from the main class
weights – weight vector [w_i,….,w_{N_s}] where w_i = area_cell_i/area_grid. Only needed if grid is non-uniform & MEMORY_SAVING== True
- Returns:
K (: np.array) if the memory saving is not active. None type otherwise.
- modulo_vki.core._pod_time.Temporal_basis_POD(K, SAVE_T_POD=False, FOLDER_OUT='./', n_Modes=10, eig_solver: str = 'eigh')[source]¶
This method computes the POD basis. For some theoretical insights, you can find the theoretical background of the proper orthogonal decomposition in a nutshell here: https://youtu.be/8fhupzhAR_M
- Parameters:
FOLDER_OUT – str. Folder in which the results will be saved (if SAVE_T_POD=True)
K – np.array. Temporal correlation matrix
SAVE_T_POD – bool. A flag deciding whether the results are saved on disk or not. If the MEMORY_SAVING feature is active, it is switched True by default.
n_Modes – int. Number of modes that will be computed
svd_solver – str. Svd solver to be used throughout the computation
- Returns:
Psi_P: np.array. POD’s Psis
- Returns:
Sigma_P: np.array. POD’s Sigmas
- modulo_vki.core._pod_space.Spatial_basis_POD(D, PSI_P, Sigma_P, MEMORY_SAVING, N_T, FOLDER_OUT='./', N_PARTITIONS=1, SAVE_SPATIAL_POD=False, rescale=False)[source]¶
This function computs the POD spatial basis from the temporal basis,
- Parameters:
D – np.array. matrix on which to project the temporal basis
PSI_P – np.array. POD’s Psis
Sigma_P – np.array. POD’s Sigmas
MEMORY_SAVING – bool. Inherited from main class, if True turns on the MEMORY_SAVING feature, loading the partitions and starting the proper algorithm
N_T – int. Number of temporal snapshots
FOLDER_OUT – str. Folder in which the results are saved if SAVE_SPATIAL_POD = True
N_PARTITIONS – int. Number of partitions to be loaded. If D has been partitioned using MODULO, this parameter is automatically inherited from the main class. To be specified otherwise.
SAVE_SPATIAL_POD – bool. If True, results are saved on disk and released from memory
rescale – bool. If False, the Sigmas are used for the normalization. If True, these are ignored and the normalization is carried out. For the standard POD, False is the way to go. However, for other decompositions (eg. the SPOD_s) you must use rescale=True
- Return Phi_P:
np.array. POD’s Phis
Functions for the mPOD¶
Similarly, the equivalent version of these functions for the mPOD are:
- modulo_vki.core._mpod_time.temporal_basis_mPOD(K, Nf, Ex, F_V, Keep, boundaries, MODE='reduced', dt=1, FOLDER_OUT: str = './', MEMORY_SAVING: bool = False, SAT: int = 100, n_Modes=10, eig_solver: str = 'svd_sklearn_randomized')[source]¶
This function computes the PSIs for the mPOD. In this implementation, a “dft-trick” is proposed, in order to avoid expansive SVDs. Randomized SVD is used by default for the diagonalization.
- Parameters:
K – np.array Temporal correlation matrix
dt – float. 1/fs, the dt between snapshots. Units in seconds.
Nf – np.array. Vector collecting the order of the FIR filters used in each scale.
Ex – int. Extension at the boundaries of K to impose the boundary conditions (see boundaries). It must be at least as Nf.
F_V – np.array. Frequency splitting vector, containing the frequencies of each scale (see article). If the time axis is in seconds, these frequencies are in Hz.
Keep – np.array. Vector defining which scale to keep.
boundaries – str -> {‘nearest’, ‘reflect’, ‘wrap’ or ‘extrap’}. In order to avoid ‘edge effects’ if the time correlation matrix is not periodic, several boundary conditions can be used. Options are (from scipy.ndimage.convolve): ‘reflect’ (d c b a | a b c d | d c b a) The input is extended by reflecting about the edge of the last pixel. ‘nearest’ (a a a a | a b c d | d d d d) The input is extended by replicating the last pixel. ‘wrap’ (a b c d | a b c d | a b c d) The input is extended by wrapping around to the opposite edge.
MODE – tr -> {‘reduced’, ‘complete’, ‘r’, ‘raw’} As a final step of this algorithm, the orthogonality is imposed via a QR-factorization. This parameterd define how to perform such factorization, according to numpy. Options: this is a wrapper to np.linalg.qr(_, mode=MODE). Check numpy’s documentation. if ‘reduced’ The final basis will not necessarely be full. If ‘complete’ The final basis will always be full
FOLDER_OUT – str. This is the directory where intermediate results will be stored if the memory saving is active.It will be ignored if MEMORY_SAVING=False.
MEMORY_SAVING – Bool. If memory saving is active, the results will be saved locally. Nevertheless, since Psi_M is usually not expensive, it will be returned.
SAT – int. Maximum number of modes per scale. The user can decide how many modes to compute; otherwise, modulo set the default SAT=100.
n_Modes – int. Total number of modes that will be finally exported
eig_solver – str. This is the eigenvalue solver that will be used. Refer to eigs_swith for the options.
- Return PSI_M:
np.array. The mPOD PSIs. Yet to be sorted !
- modulo_vki.core._mpod_space.spatial_basis_mPOD(D, PSI_M, N_T, N_PARTITIONS, N_S, MEMORY_SAVING, FOLDER_OUT, SAVE: bool = False, weights: array = array([], dtype=float64))[source]¶
Given the temporal basis of the mPOD now the spatial ones are computed
- Parameters:
D – Snapshot matrix D: if memory savig is active, this is ignored.
PSI_M – np.array.: The mPOD temporal basis Psi tentatively assembled from all scales
N_T – int. Number of snapshots
N_PARTITIONS – int. Number of partitions in the memory saving
N_S – int. Number of grid points in space
MEMORY_SAVING – bool. Inherited from main class, if True turns on the MEMORY_SAVING feature, loading the partitions and starting the proper algorithm
FOLDER_OUT – str. Folder in which the results are saved if SAVE_SPATIAL_POD = True
SAVE_SPATIAL_POD – bool. If True, results are saved on disk and released from memory
weights – np.array weight vector [w_i,….,w_{N_s}] where w_i = area_cell_i/area_grid. Only needed if grid is non-uniform & MEMORY_SAVING== True
- Returns:
Phi_M, Psi_M, Sigma_M: np.arrays. The final (sorted) mPOD decomposition
Functions for the DFT¶
- modulo_vki.core._dft.dft_fit(N_T, F_S, D, FOLDER_OUT, SAVE_DFT=False)[source]¶
This function computes the DFT form the dataset D. Currently, this does not handle the memory saving feature.
- Parameters:
N_T – int. number of snapshots
F_S – Sampling frequency (in Hz)
D – Snapshot matrix
FOLDER_OUT – Folder in which the results are saved if SAVE_SPATIAL_POD = True
SAVE_DFT – If True, results are saved on disk and released from memory
- Returns:
Sorted_Freqs, np.array Frequency bins, in Hz.
- Returns:
Phi_F, np.array (Complex) Spatial structures for each mode
- Returns:
SIGMA_F, np.array (real) amplitude of each modes
Functions for the DMD (PIP)¶
- modulo_vki.core._dmd_s.dmd_s(D_1, D_2, n_Modes, F_S, SAVE_T_DMD=False, FOLDER_OUT='./', svd_solver: str = 'svd_sklearn_truncated')[source]¶
This method computes the Dynamic Mode Decomposition (DMD) using hte PIP algorithm from Penland.
- Parameters:
D_1 – np.array First portion of the data, i.e. D[:,0:n_t-1]
D_2 – np.array Second portion of the data, i.e. D[:,1:n_t]
Sigma_P (Phi_P, Psi_P,) – np.arrays POD decomposition of D1
F_S – float Sampling frequency in Hz
FOLDER_OUT – str Folder in which the results will be saved (if SAVE_T_DMD=True)
K – np.array Temporal correlation matrix
SAVE_T_POD – bool A flag deciding whether the results are saved on disk or not. If the MEMORY_SAVING feature is active, it is switched True by default.
n_Modes – int number of modes that will be computed
svd_solver – str, svd solver to be used
- Return1 Phi_D:
np.array. DMD’s complex spatial structures
- Return2 Lambda_D:
np.array. DMD Eigenvalues (of the reduced propagator)
- Return3 freqs:
np.array. Frequencies (in Hz, associated to the DMD modes)
- Return4 a0s:
np.array. Initial Coefficients of the Modes
For the SPODs, currently these are implemented in a single function. MODULO acts as a wrapper to these. The memory saving feature on these is not yet implemented.
Functions for the SPOD_S¶
- modulo_vki.core._spod_s.compute_SPOD_s(D, K, F_S, n_s, n_t, N_o=100, f_c=0.3, n_Modes=10, SAVE_SPOD=True, FOLDER_OUT='./', MEMORY_SAVING=False, N_PARTITIONS=1)[source]¶
This method computes the Spectral POD of your data. This is the one by Sieber et al (https://www.cambridge.org/core/journals/journal-of-fluid-mechanics/article/abs/spectral-proper-orthogonal-decomposition/DCD8A6EDEFD56F5A9715DBAD38BD461A)
- Parameters:
F_S – float, Sampling Frequency [Hz]
N_o – float, Semi-Order of the diagonal filter. Note that the filter order will be 2 N_o +1 (to make sure it is odd)
f_c – float, cut-off frequency of the diagonal filter
n_Modes – float, number of modes to be computed
SAVE_SPOD – bool, If True, MODULO will save the output in self.FOLDER OUT/MODULO_tmp
FOLDER_OUT – string Define where the out will be stored (ignored if SAVE_POD=False)
SAVING (MEMORY) – bool Define if memory saving is active or not (reduntant; to be improved) Currently left for compatibility with the rest of MODULO.
N_PARTITIONS – int number of partitions (if memory saving = False, it should be 1). (reduntant; to be improved) Currently left for compatibility with the rest of MODULO.
- Return Psi_P:
np.array SPOD Psis
- Return Sigma_P:
np.array SPOD Sigmas.
- Return Phi_P:
np.array SPOD Phis
Functions for the SPOD_T¶
- modulo_vki.core._spod_t.compute_SPOD_t(D, F_S, L_B=500, O_B=250, n_Modes=10, SAVE_SPOD=True, FOLDER_OUT='/', possible_svds='svd_sklearn_truncated')[source]¶
This method computes the Spectral POD of your data. This is the one by Town et al (https://www.cambridge.org/core/journals/journal-of-fluid-mechanics/article/spectral-proper-orthogonal-decomposition-and-its-relationship-to-dynamic-mode-decomposition-and-resolvent-analysis/EC2A6DF76490A0B9EB208CC2CA037717)
- Parameters:
D – array. snapshot matrix to decompose, of size N_S,N_T
F_S – float, Sampling Frequency [Hz]
L_B – float, Lenght of the chunks
O_B – float, Overlapping between blocks in the chunk
n_Modes – float, Number of modes to be computed FOR EACH FREQUENCY
SAVE_SPOD – bool, If True, MODULO will save the output in FOLDER OUT/MODULO_tmp
possible_svds – str, Svd solver to be used throughout the computation
- Return Psi_P_hat:
np.array Spectra of the SPOD Modes
- Return Sigma_P:
np.array Amplitudes of the SPOD Modes.
- Return Phi_P:
np.array SPOD Phis
- Return freq:
float Frequency bins for the Spectral POD