pyserep.analysis

pyserep.analysis — validation, performance, convergence, sensitivity.

class pyserep.analysis.ConvergencePoint(param_value: float, n_modes: int, n_dofs: int, kappa: float, max_frf_err_pct: float, rms_frf_err_pct: float, max_freq_err_pct: float)[source]

Bases: object

Single point in a convergence sweep.

kappa: float
max_freq_err_pct: float
max_frf_err_pct: float
n_dofs: int
n_modes: int
param_value: float
rms_frf_err_pct: float
class pyserep.analysis.ConvergenceStudy(param_name: str, param_label: str, points: List[ConvergencePoint])[source]

Bases: object

Results from a full convergence sweep.

param_label: str
param_name: str
plot(save_path: str | None = None, show: bool = False) None[source]

Plot FRF error and condition number vs parameter.

points: List[ConvergencePoint]
table() str[source]

Return the convergence results as a formatted ASCII table string.

class pyserep.analysis.PerformanceMetrics(n_full_dofs: int, n_selected_modes: int, n_master_dofs: int, dof_reduction_pct: float, mode_retention_pct: float, kappa: float, frf_method: str, n_freq_points: int, n_bands: int, frf_flops_rom: int, frf_flops_ref: int, frf_speedup: float, t_eigensolver_s: float, t_mode_select_s: float, t_dof_select_s: float, t_rom_build_s: float, t_frf_s: float, t_total_s: float)[source]

Bases: object

All performance and efficiency metrics for a pipeline run.

dof_reduction_pct: float
frf_flops_ref: int
frf_flops_rom: int
frf_method: str
frf_speedup: float
kappa: float
mode_retention_pct: float
n_bands: int
n_freq_points: int
n_full_dofs: int
n_master_dofs: int
n_selected_modes: int
summary() str[source]

Return a formatted multi-line performance summary string.

t_dof_select_s: float
t_eigensolver_s: float
t_frf_s: float
t_mode_select_s: float
t_rom_build_s: float
t_total_s: float
class pyserep.analysis.ValidationReport(eigenvalue_errors_pct: ndarray, max_eigenvalue_error_pct: float, mean_eigenvalue_error_pct: float, mass_ortho_error: float, stiff_ortho_error: float, expansion_error: float, mac_values: ndarray, min_mac: float, mean_mac: float, ka_positive_definite: bool, ma_positive_definite: bool, ka_condition_number: float, ma_condition_number: float)[source]

Bases: object

Full validation report for a SEREP ROM.

eigenvalue_errors_pct: ndarray
expansion_error: float
ka_condition_number: float
ka_positive_definite: bool
ma_condition_number: float
ma_positive_definite: bool
mac_values: ndarray
mass_ortho_error: float
max_eigenvalue_error_pct: float
mean_eigenvalue_error_pct: float
mean_mac: float
min_mac: float
passed(strict: bool = False) bool[source]

Return True if all key checks pass.

stiff_ortho_error: float
summary() str[source]

Return a formatted multi-line validation report as a string.

pyserep.analysis.dof_count_study(K: csc_matrix, M: csc_matrix, phi: ndarray, freqs_hz: ndarray, selected_modes: ndarray, force_dofs: List[int], output_dofs: List[int], n_master_values: List[int], freq_eval: ndarray, H_ref: dict, zeta: float = 0.001, verbose: bool = True) ConvergenceStudy[source]

Study how FRF accuracy and condition number change with n_master.

Fixes the mode set and varies the number of master DOFs from len(selected_modes) up to max(n_master_values).

Parameters:
  • K (full model)

  • M (full model)

  • phi (full model)

  • freqs_hz (full model)

  • selected_modes (np.ndarray of int — fixed mode set)

  • force_dofs (list of int)

  • output_dofs (list of int)

  • n_master_values (list of int — DOF counts to sweep)

  • freq_eval (np.ndarray — frequency grid for FRF)

  • H_ref (dict — reference FRF on freq_eval)

  • zeta (float)

  • verbose (bool)

Return type:

ConvergenceStudy

pyserep.analysis.eigenvalue_error(Ka: ndarray, Ma: ndarray, full_freqs_hz: ndarray, selected_modes: ndarray, verbose: bool = True) Tuple[ndarray, float][source]

Compute per-mode eigenvalue preservation error (%).

Returns:

  • errors_pct (np.ndarray, shape (m,))

  • max_error_pct (float)

pyserep.analysis.eigenvalue_sensitivity(K: spmatrix, M: spmatrix, phi: ndarray, freqs_hz: ndarray, selected_modes: ndarray, dK_dp: spmatrix, dM_dp: spmatrix, verbose: bool = True) ndarray[source]

Compute eigenvalue sensitivities ∂λᵢ/∂p using Nelson’s method.

For mass-normalised modes (φᵢᵀ M φᵢ = 1):

∂λᵢ/∂p = φᵢᵀ (∂K/∂p − λᵢ ∂M/∂p) φᵢ

Parameters:
  • K (sparse matrices — full structural matrices)

  • M (sparse matrices — full structural matrices)

  • phi (np.ndarray, shape (N, n_all_modes))

  • freqs_hz (np.ndarray — full-model natural frequencies)

  • selected_modes (np.ndarray of int)

  • dK_dp (sparse matrices — parameter derivative matrices) – Represents ∂K/∂p and ∂M/∂p for a single parameter p. For a Young’s modulus perturbation: dK_dp = K/E, dM_dp = 0.

  • dM_dp (sparse matrices — parameter derivative matrices) – Represents ∂K/∂p and ∂M/∂p for a single parameter p. For a Young’s modulus perturbation: dK_dp = K/E, dM_dp = 0.

  • verbose (bool)

Returns:

dlam_dp

Return type:

np.ndarray, shape (m,) — ∂λᵢ/∂p in (rad/s)² per unit of p

Examples

Sensitivity of eigenvalues to a 1% stiffness increase:

>>> dK_dp = 0.01 * K   # 1% increase in E
>>> dM_dp = sp.csc_matrix(K.shape)   # mass unchanged
>>> dlam = eigenvalue_sensitivity(K, M, phi, freqs, modes, dK_dp, dM_dp)
>>> dfreq_hz = dlam / (2 * np.pi * freqs_hz[modes]) / (2 * np.pi)
pyserep.analysis.flop_count(n_modes: int, n_freq: int, n_pairs: int, method: str = 'direct') int[source]

Estimate FRF computation FLOP count.

Parameters:
  • n_modes (int) – Number of retained modes (m for ROM, all elastic for reference).

  • n_freq (int) – Number of evaluation frequency points.

  • n_pairs (int) – Number of force/output DOF pairs.

  • method (str) – "direct" or "modal".

Returns:

Approximate floating-point operation count.

Return type:

int

Notes

Modal FRF: per frequency per mode per pair: ~8 FLOPs (2 multiplications + complex division + addition)

Direct FRF: per frequency: one m×m complex LU factor + m×m backsolve LU factor: (2/3) m³ FLOPs Backsolve: m² FLOPs per right-hand side × n_pairs

pyserep.analysis.frf_sensitivity(Ka: ndarray, Ma: ndarray, dKa_dp: ndarray, dMa_dp: ndarray, local_force: List[int], local_output: List[int], freq_eval: ndarray, zeta: float = 0.001, damping_type: str = 'modal') Dict[str, ndarray][source]

Compute the FRF parameter sensitivity ∂H/∂p via direct differentiation.

Differentiating Z(ω) H = I with respect to p gives:

∂H/∂p = −Z⁻¹ (∂Z/∂p) Z⁻¹

where ∂Z/∂p = ∂Kₐ/∂p − ω² ∂Mₐ/∂p (for viscous damping proportional to Ka, Ma).

Parameters:
  • Ka (np.ndarray, shape (m, m) — reduced matrices)

  • Ma (np.ndarray, shape (m, m) — reduced matrices)

  • dKa_dp (np.ndarray, shape (m, m) — parameter derivatives of Ka, Ma)

  • dMa_dp (np.ndarray, shape (m, m) — parameter derivatives of Ka, Ma)

  • local_force (list of int — local DOF indices within master set)

  • local_output (list of int — local DOF indices within master set)

  • freq_eval (np.ndarray — evaluation frequencies (Hz))

  • zeta (damping parameters (same as compute_frf_direct))

  • damping_type (damping parameters (same as compute_frf_direct))

Returns:

dH_dp – Values are complex arrays of shape (n_freq,) representing ∂H/∂p.

Return type:

dict — same key structure as compute_frf_direct()

Notes

This gives the first-order sensitivity. For a parameter change Δp, the FRF changes approximately by ΔH ≈ (∂H/∂p) Δp.

pyserep.analysis.material_perturbation_study(Ka_nominal: ndarray, Ma_nominal: ndarray, local_force: List[int], local_output: List[int], freq_eval: ndarray, H_nominal: Dict[str, ndarray], param_values: List[float], Ka_func: Callable[[float], ndarray], Ma_func: Callable[[float], ndarray], zeta: float = 0.001, verbose: bool = True) Dict[str, ndarray][source]

Compute FRFs for a sweep of material parameter values.

Useful for comparing Standard vs HyAPA material configurations, or for uncertainty quantification (±5% E, ±3% ρ, etc.).

Parameters:
  • Ka_nominal (np.ndarray — baseline reduced matrices)

  • Ma_nominal (np.ndarray — baseline reduced matrices)

  • local_force (list of int)

  • local_output (list of int)

  • freq_eval (np.ndarray)

  • H_nominal (dict — baseline FRF (from compute_frf_direct))

  • param_values (list of float — parameter values to evaluate) – e.g. [0.95, 0.97, 1.00, 1.03, 1.05] for ±5% sweep

  • Ka_func (callable — Ka_func(p) → np.ndarray(m, m)) – Returns Ka for parameter value p. For stiffness scaling: lambda p: p * Ka_nominal

  • Ma_func (callable — Ma_func(p) → np.ndarray(m, m)) – For mass scaling: lambda p: p * Ma_nominal

  • zeta (see compute_frf_direct)

  • verbose (see compute_frf_direct)

Returns:

"param_values" — np.ndarray of the swept values "H_sweep" — np.ndarray, shape (n_params, n_freq) — FRF magnitudes "max_deviation_pct" — np.ndarray, shape (n_params,)

Return type:

dict with keys

pyserep.analysis.modal_assurance_criterion(phi_ref: ndarray, phi_rom: ndarray) ndarray[source]

Compute the full MAC matrix between two mode shape sets.

MAC[i,j] = |φᵢᵀ φⱼ|² / (|φᵢ|² |φⱼ|²)

Parameters:
  • phi_ref (np.ndarray, shape (N, m))

  • phi_rom (np.ndarray, shape (N, m))

Returns:

mac_matrix – Values in [0, 1]. Diagonal near 1 → good mode pairing.

Return type:

np.ndarray, shape (m, m)

pyserep.analysis.mode_count_study(K: csc_matrix, M: csc_matrix, phi: ndarray, freqs_hz: ndarray, force_dofs: List[int], output_dofs: List[int], f_max: float, f_max_values: List[float], zeta: float = 0.001, rb_hz: float = 1.0, n_freq: int = 500, verbose: bool = True) ConvergenceStudy[source]

Study FRF convergence as the upper frequency cutoff is varied.

For each value in f_max_values, runs a full SEREP pipeline (mode selection → DOF selection → ROM build → direct FRF) and records accuracy metrics.

Parameters:
  • K (sparse matrices)

  • M (sparse matrices)

  • phi (modal matrix and frequencies (pre-computed))

  • freqs_hz (modal matrix and frequencies (pre-computed))

  • force_dofs (list of int)

  • output_dofs (list of int)

  • f_max (float) – Upper limit of the reference FRF (Hz).

  • f_max_values (list of float) – Cutoff frequencies to sweep (Hz). Each defines the upper edge of a single analysis band.

  • zeta (float)

  • rb_hz (float)

  • n_freq (int) – Number of FRF evaluation points.

  • verbose (bool)

Return type:

ConvergenceStudy

pyserep.analysis.monte_carlo_frf(Ka_nominal: ndarray, Ma_nominal: ndarray, local_force: List[int], local_output: List[int], freq_eval: ndarray, E_cov_pct: float = 2.0, rho_cov_pct: float = 1.0, n_samples: int = 50, seed: int = 42, zeta: float = 0.001, verbose: bool = True) Dict[str, ndarray][source]

Monte Carlo FRF uncertainty quantification.

Perturbs Ka (proportional to E uncertainty) and Ma (proportional to ρ uncertainty) and computes the FRF for each sample, returning the mean, standard deviation, and 5th/95th percentile bands.

Parameters:
  • Ka_nominal (np.ndarray, shape (m, m))

  • Ma_nominal (np.ndarray, shape (m, m))

  • local_force (list of int)

  • local_output (list of int)

  • freq_eval (np.ndarray)

  • E_cov_pct (float — coefficient of variation for Young's modulus (%))

  • rho_cov_pct (float — coefficient of variation for density (%))

  • n_samples (int)

  • seed (int)

  • zeta

  • verbose

Returns:

"H_mean" — mean FRF magnitude "H_std" — standard deviation "H_p5" — 5th percentile "H_p95" — 95th percentile "H_all" — all samples, shape (n_samples, n_freq)

Return type:

dict with keys

pyserep.analysis.orthogonality_check(phi: ndarray, M: csc_matrix, selected_modes: ndarray) float[source]

Return |ΦᵀMΦ I|_max. Should be < 1e-10 for mass-normalised modes.

pyserep.analysis.reduction_metrics(n_full_dofs: int, n_master_dofs: int, n_all_modes: int, n_selected_modes: int) Dict[str, float][source]

Compute reduction ratio metrics.

Returns:

dof_reduction_ratio, dof_retention_pct, mode_retention_pct, size_ratio

Return type:

dict with keys

pyserep.analysis.summarise_performance(n_full_dofs: int, n_selected_modes: int, n_master_dofs: int, n_all_modes: int, kappa: float, n_freq: int, n_bands: int, n_pairs: int, frf_method: str = 'direct', frf_flops_rom: int | None = None, frf_flops_ref: int | None = None, t_eigensolver_s: float = 0.0, t_mode_select_s: float = 0.0, t_dof_select_s: float = 0.0, t_rom_build_s: float = 0.0, t_frf_s: float = 0.0, t_total_s: float = 0.0) PerformanceMetrics[source]

Collect all performance metrics into a PerformanceMetrics.

Parameters are self-explanatory from the return type documentation.

pyserep.analysis.validate_serep(K: csc_matrix, M: csc_matrix, phi: ndarray, freqs_hz: ndarray, selected_modes: ndarray, master_dofs: ndarray, T: ndarray, Ka: ndarray, Ma: ndarray, verbose: bool = True) ValidationReport[source]

Run the complete SEREP validation suite.

Parameters:
  • K (sparse matrices)

  • M (sparse matrices)

  • phi (np.ndarray, shape (N, n_modes))

  • freqs_hz (np.ndarray)

  • selected_modes (np.ndarray of int)

  • master_dofs (np.ndarray of int)

  • T (np.ndarray, shape (N, m) — transformation matrix)

  • Ka (np.ndarray, shape (m, m) — reduced matrices)

  • Ma (np.ndarray, shape (m, m) — reduced matrices)

  • verbose (bool)

Return type:

ValidationReport