MATLAB/Octave Driver#
The main driver script is drive/matlab/rom_online_solver.m, which loads the basis functions and runs the ROM in MATLAB/Octave using a BDF3/EXT3 time-stepper. Supporting functions are defined in separate files and documented below.
Plotting#
The Matlab/Octave driver uses NekToolKit for plotting of 2D fields. To enable this, clone the NekToolKit repository and either append it to the MATLABPATH environment variable (for MATLAB) or call addpath within octaverc (for Octave).
git clone https://github.com/kent0/NekToolKit
export MATLABPATH=$(pwd)/NekToolKit/matlab
Operators#
- conv_deim(ucoef, rom_data, method)#
APPLY_ROM_CONVECTION Unified runtime for DEIM, CLSDEIM, and MCLSDEIM method: ‘deim’, ‘cls’, or ‘mcls’
- conv_fom(ucoef, pod_u, pod_v, x, y, varargin)#
Optimized Hybrid Pseudo-ROM / FOM Convection
- conv_tensor_dense(ucoef, pod_u, pod_v, x, y, tensor_size)#
C_ijk = < phi_k , (phi_i . grad) phi_j >
- conv_tensor(ucoef, pod_u, pod_v, x, y)#
Optimized Convection Tensor for NekROM C(i,j,k) = < phi_i . grad(phi_j) , phi_k >
- conv_tensor_sparse(ucoef, pod_u, pod_v, x, y, validator)#
De-aliased Sparse Convection Tensor with Automated Symmetric Validation
- gen_Au(pod_u, pod_v, x, y)#
GEN_AU SEM-consistent ROM diffusion and mass operators Compatible with MATLAB and GNU Octave
- get_Me(x, y)#
- lgrad(u, x, y, mode, reset)#
LGRAD Computes the local gradient on a spectral element. Optional ‘reset’ argument forces recalculation of geometry.
- lcurl(u, v, x, y)#
Input and Output#
- get_grid(snaps, reorder)#
Handle default value for reorder
- get_pod_basis_from_arrays(u_snaps, v_snaps, x, y, nb, subtract_mean, conserve_momentum, method, inner_product)#
GET_POD_BASIS_FROM_ARRAYS Performs POD with mass-matrix or diffusion weighting Updated for Octave/MATLAB compatibility (2026)
- get_pod_basis(snaps_obj, nb, reorder, subtract_mean, conserve_momentum, method, inner_product)#
Assemble the snapshots
- get_r_dim_ops(au_full, bu_full, cu_full, u0_full, uk_full, nb)#
Linear Dynamics extraction
- get_snaps(snaps, reorder)#
Default reorder to 0 if not provided
- get_sort_order(x, y)#
1. Calculate means by nesting (compatible with Octave & MATLAB) We mean across dim 1, then dim 2
- load_full_ops(path)#
load_full_ops function loads the full ROM operators stored in the specified path specified by input variable path.
Output: a0_full : Full stiffness matrix of size mb+1 x mb+1 (The +1 comes from the zeroth mode) b0_full : Full mass matrix of size mb+1 x mb+1 cu_full : Full advection tensor of size mb x mb+1 x mb+1 u0_full : Vector of size mb+1, which contains the ROM coefficients of the
projection of initial conditons onto mb-reduced space
- uk_fullMatrix of size mb+1 x ns. Each column contains the ROM coefficients
of the projection of the snapshot onto mb-reduced space
mb : total number of modes ns : nubmer of snapshots used to create those operators
This function can take times to load operators with nb >= 300.
- output_fields(basepath, data, ifvort, ifwrite, ifvis)#
Pre-fetch size to avoid repeated struct access
- write_field(basepath, data, wdsize)#
Determined Dimensions
DEIM Point Generators#
- s_opt(Vo, M, index, outfile)#
S_OPT Generates S-optimal indices of a matrix Optimized implementation focusing on vectorization and memory efficiency. Obtained by querying Gemini for improvements to the original algorithm. Original code in the extra folder.
- Parameters:
Vo – Candidate Matrix (N_Bm x N)
M – Number of S-optimal points to compute
index – Initial index set (use [] if none)
outfile – Optional filename string for saving results
- gpode()#
The type of the None singleton.
- gnat(Q, m_used, nsr)#
- gappy_pod(U_nl, n, method, tol)#
GAPPY_POD_QR Optimized Gappy POD / DEIM with oversampling.
- USAGE:
[indices] = gappy_pod(U_nl, n, ‘synced’, 1e-10)
- INPUTS:
U_nl - Basis matrix (N x p), usually POD modes. n - Maximum number of indices to select. method - ‘synced’ (default) or ‘spread’. tol - Early exit threshold for residual norm.