Skip to main content
Ctrl+K

femorph-solver

  • Getting started
  • User guide
  • Reference
  • Examples
  • API reference
    • Changelog
    • Provenance inventory
  • GitHub
  • Getting started
  • User guide
  • Reference
  • Examples
  • API reference
  • Changelog
  • Provenance inventory
  • GitHub

Section Navigation

  • FEA in three phases
  • Pre-processing
    • The Model
    • Element library
      • TRUSS2 — 3-D 2-node axial bar
      • SPRING — 3-D 2-node longitudinal spring
      • POINT_MASS — concentrated nodal mass
      • BEAM2 — 3D 2-node Euler–Bernoulli beam
      • QUAD4_PLANE — 2-D 4-node structural solid
      • QUAD4_SHELL — 4-node Mindlin–Reissner shell
      • HEX8 — 8-node hexahedral solid
      • HEX20 — 20-node quadratic hexahedral solid
      • TET10 — 10-node quadratic tetrahedron
    • Materials
    • Real constants and sections
    • Boundary conditions and loads
  • Solving
    • Static analysis
    • Modal analysis
    • Transient analysis
    • Cyclic-symmetry modal analysis
    • Solvers
    • Performance
    • Thread control
    • Benchmarking femorph-solver on your machine
    • Estimating wall time + memory before you solve
    • Solver exploration (extras)
    • Out-of-core (OOC) Pardiso — limits and tuning
  • Post-processing
    • Result objects
    • Strain and stress recovery
    • Visualisation and export
  • MAPDL compatibility
  • Cross-solver comparisons
  • User guide
  • Solving
  • Modal analysis

Modal analysis#

femorph_solver.solvers.modal.solve_modal() returns the lowest free-vibration modes of a structure. Mathematically, it solves the generalised eigenproblem

\[K \, \varphi = \omega^2 \, M \, \varphi\]

for the smallest \(n_\text{modes}\) eigenvalues \(\omega^2\) and their eigenvectors \(\varphi\) on the free-DOF partition of the model.

Dense vs sparse path#

Two paths, chosen by the number of free DOFs:

  • Dense LAPACK (scipy.linalg.eigh) — used automatically when n_free <= DENSE_CUTOFF (currently 500). Fastest on tiny problems; returns every mode.

  • Sparse eigen backend — used otherwise; selected via eigen_solver="arpack" | "lobpcg" | "primme" | "dense". Default is "arpack" with shift-invert.

The shift-invert step inside ARPACK reuses the linear-solver registry from Solvers — pass linear_solver="cholmod" (or any other name) to steer the inner factorisation.

API#

from femorph_solver.solvers.modal import solve_modal

K = m.stiffness_matrix()
M = m.mass_matrix()

# Clamped at node 1: the pin produces 3 Dirichlet DOFs.
prescribed = {int(i): 0.0 for i in clamped_dofs}

result = solve_modal(
    K, M,
    prescribed=prescribed,
    n_modes=10,
    eigen_solver="arpack",
    linear_solver="auto",
    sigma=0.0,
    tol=0.0,
)

result.omega_sq                   # (n_modes,) ω²
result.frequency                  # (n_modes,) f [Hz] = √ω² / 2π
result.mode_shapes                # (N, n_modes), full-length vectors
result.free_mask                  # (N,) bool — DOFs kept in the eigenproblem

Or m.modal_solve(n_modes=10) as a shortcut that pulls prescribed from the Model’s D records.

The shift parameter sigma#

ARPACK’s shift-invert path solves \((K - \sigma M)^{-1} M \varphi = (1/(\omega^2 - \sigma)) \varphi\) and returns eigenvalues nearest \(\sigma\). Useful in two cases:

  1. Free-free structures (no Dirichlet). \(K\) has a six-dimensional rigid-body kernel, and sigma=0 would try to factor a singular matrix. Pass sigma=-1.0 (or any small negative number) so \(K - \sigma M\) is SPD:

    solve_modal(K, M, prescribed={}, n_modes=10,
                eigen_solver="arpack", sigma=-1.0)
    
  2. Targeting a frequency band. If you want modes near 2 kHz, pass sigma = (2 * π * 2000)² ≈ 1.58e8. ARPACK converges toward those eigenvalues instead of the smallest absolute ones.

Tolerance#

tol=0.0 asks for machine precision (the ARPACK default). Raise the value only when you need very-fast-but-imprecise iteration — the relative error on an eigenvalue scales roughly linearly with tol.

Thread control#

The solver accepts thread_limit: int | None = None and wraps the entire dense / sparse eigenpath in the thread-cap context. See Thread control.

Rigid-body DOFs#

Just like Static analysis, any DOF with |K_ii| <= 1e-12 * max(|K_ii|) is treated as constrained and dropped from the eigenproblem. That lets a pure-axial truss’s transverse DOFs fall out automatically without the user having to find and fix them.

See also#

  • Cyclic-symmetry modal analysis — cyclic-symmetry modal (sector-only reduction).

  • Solvers — sparse linear backend for the shift-invert inner solve.

  • Result objects — reading the ModalResult.

On this page
  • Dense vs sparse path
  • API
  • The shift parameter sigma
  • Tolerance
  • Thread control
  • Rigid-body DOFs
  • See also
Show Source

© Copyright 2026, Alex Kaszynski.

Created using Sphinx 9.1.0.

Built with the PyData Sphinx Theme 0.17.1.