pyserep.core

pyserep.core — eigensolver and SEREP ROM builder.

pyserep.core.build_rayleigh_damping(Ka: ndarray, Ma: ndarray, zeta: float, freqs_hz: ndarray, modes: ndarray) ndarray[source]

Build a Rayleigh damping matrix for the reduced model.

Rayleigh damping: Cₐ = α Mₐ + β Kₐ

where α and β are chosen to give damping ratio zeta at the first and last selected natural frequencies.

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

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

  • zeta (float) – Desired damping ratio (uniform across modes).

  • freqs_hz (np.ndarray) – Full-model natural frequencies.

  • modes (np.ndarray of int) – Selected mode indices.

Returns:

Ca – Rayleigh damping matrix.

Return type:

np.ndarray, shape (m, m)

pyserep.core.build_serep_rom(K: csc_matrix, M: csc_matrix, phi: ndarray, selected_modes: ndarray, master_dofs: ndarray, verbose: bool = True) Tuple[ndarray, ndarray, ndarray][source]

Build the SEREP transformation matrix T and compute reduced matrices Kₐ, Mₐ.

Parameters:
  • K (scipy.sparse.csc_matrix) – Full stiffness matrix (N × N). Must be real and symmetric.

  • M (scipy.sparse.csc_matrix) – Full mass matrix (N × N). Must be real, symmetric, and positive definite.

  • phi (np.ndarray, shape (N, n_all_modes)) – Full mass-normalised modal matrix.

  • selected_modes (np.ndarray of int) – Mode indices to retain (output of select_modes()).

  • master_dofs (np.ndarray of int) – Master DOF indices (output of select_dofs_eid()). Must satisfy len(master_dofs) == len(selected_modes) for exact SEREP.

  • verbose (bool) – Print construction diagnostics.

Returns:

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

  • Ka (np.ndarray, shape (m, m)) – Reduced stiffness matrix.

  • Ma (np.ndarray, shape (m, m)) – Reduced mass matrix.

Raises:
  • ValueError – If len(master_dofs) != len(selected_modes) and the exact inverse would be undefined.

  • RuntimeWarning – If the condition number of Φₐ is very large (> 10⁶), indicating a poorly conditioned master DOF set.

Notes

When len(master_dofs) > len(selected_modes) (over-constrained), the Moore–Penrose pseudoinverse is used, resulting in a least-squares SEREP approximation.

Examples

>>> T, Ka, Ma = build_serep_rom(K, M, phi, selected_modes, master_dofs)
>>> Ka.shape
(37, 37)
pyserep.core.solve_eigenproblem(K: csc_matrix, M: csc_matrix, n_modes: int, sigma: float = 0.01, tol: float = 1e-10, maxiter: int | None = None, verbose: bool = True) Tuple[ndarray, ndarray][source]

Solve the generalised eigenvalue problem K φ = λ M φ for the n_modes lowest eigenpairs.

Parameters:
  • K (scipy.sparse.csc_matrix) – Stiffness matrix (N × N), real symmetric positive semi-definite. Asymmetric K will produce complex eigenvalues and incorrect results.

  • M (scipy.sparse.csc_matrix) – Mass matrix (N × N), real symmetric positive definite. M must be non-singular; singular M causes ARPACK to fail.

  • n_modes (int) – Number of eigenpairs to extract. Must be < N − 1.

  • sigma (float) – Shift parameter (rad²/s²) for shift-invert. A small positive value (e.g. 0.01) ensures rigid-body modes (λ ≈ 0) are captured.

  • tol (float) – ARPACK convergence tolerance. 1e-10 is appropriate for structural analysis; decrease to 1e-12 for very high accuracy requirements.

  • maxiter (int, optional) – Maximum ARPACK iterations. Default: 10 * N.

  • verbose (bool) – Print eigenvalue summary.

Returns:

  • freqs_hz (np.ndarray, shape (n_modes,)) – Natural frequencies in Hz, sorted ascending.

  • phi (np.ndarray, shape (N, n_modes)) – Mass-normalised mode shapes. phi[:, i] corresponds to freqs_hz[i].

Raises:

Examples

>>> freqs, phi = solve_eigenproblem(K, M, n_modes=100)
>>> freqs.shape
(100,)
>>> phi.shape
(66525, 100)
pyserep.core.verify_eigenvalues(Ka: ndarray, Ma: ndarray, full_freqs_hz: ndarray, selected_modes: ndarray, verbose: bool = True) Tuple[ndarray, float][source]

Verify SEREP’s defining property: exact eigenvalue preservation.

Solves the reduced eigenvalue problem Kₐ φ = λ Mₐ φ and compares the resulting frequencies to the full-model values.

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

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

  • full_freqs_hz (np.ndarray) – All full-model natural frequencies (Hz).

  • selected_modes (np.ndarray of int) – Indices identifying the target frequencies.

  • verbose (bool)

Returns:

  • freq_errors_pct (np.ndarray, shape (m,)) – Per-mode percentage error |f_rom - f_full| / f_full × 100.

  • max_error_pct (float) – Maximum absolute percentage error. Should be < 0.001% for a well-conditioned SEREP.

Examples

>>> errors, max_err = verify_eigenvalues(Ka, Ma, freqs_hz, selected_modes)
>>> max_err < 0.001
True