Static analysis =============== :func:`femorph_solver.solvers.static.solve_static` solves the standard linear-elastic problem .. math:: K \, u = F with Dirichlet constraints imposed at a set of DOFs. It's the fastest path in femorph-solver — one sparse factorisation, one back-solve — and the starting point for more involved analyses (static + gravity, static + thermal preload, static pre-stress for modal). Partitioning ------------ Given a prescribed set :math:`c` of DOFs with fixed values :math:`u_c` and the complementary free set :math:`f`, the system partitions as .. math:: \begin{bmatrix} K_{ff} & K_{fc} \\ K_{cf} & K_{cc} \end{bmatrix} \begin{bmatrix} u_f \\ u_c \end{bmatrix} = \begin{bmatrix} F_f \\ R_c \end{bmatrix} where :math:`R_c` is the reaction force at constrained DOFs. The solve factors :math:`K_{ff}` once, solves for .. math:: u_f = K_{ff}^{-1} \bigl(F_f - K_{fc} u_c\bigr), and then back-substitutes .. math:: R_c = K_{cf} u_f + K_{cc} u_c - F_c to recover the reactions. API --- .. code-block:: python from femorph_solver.solvers.static import solve_static K = m.stiffness_matrix() F = np.zeros(K.shape[0]) # full-length load vector F[load_dof] = -1000.0 # 1 kN somewhere prescribed = {d: 0.0 for d in fixed_dofs} # e.g. a clamped face result = solve_static(K, F, prescribed=prescribed) result.displacement # (N,) — every DOF, zeros where clamped result.reaction # (N,) — non-zero only at constrained DOFs result.free_mask # (N,) bool — True where we solved Or from a Model that already carries ``D`` and ``F`` records:: m.d(1, "UX", 0.0) # ... pin some DOFs ... m.f(10, "FZ", -1000.0) result = m.solve() # wraps solve_static with the record data Backend selection ----------------- ``solve_static`` accepts ``linear_solver="auto"`` (default) which picks the fastest available SPD direct solver in this order: PARDISO → CHOLMOD → UMFPACK → SuperLU. Override with any name from :func:`femorph_solver.solvers.linear.list_linear_solvers`:: result = solve_static(K, F, prescribed, linear_solver="cholmod") See :doc:`choosing-a-solver` for the full backend table, and :doc:`performance` for measured wall times. Thread control -------------- ``thread_limit: int | None = None`` caps BLAS / OpenMP threads for the factorisation and back-solve. ``None`` defers to the process-wide default set by :func:`femorph_solver.set_thread_limit` or the ``FEMORPH_SOLVER_NUM_THREADS`` environment variable. See :doc:`threads`. Non-zero Dirichlet ------------------ Any non-zero ``u_c`` is fine — the RHS correction :math:`F_f - K_{fc} u_c` flows through automatically, and the reaction recovery produces the force needed to hold the prescribed displacement. Zero-diagonal DOFs (rigid bodies) --------------------------------- Some element formulations produce DOFs with zero diagonal stiffness (e.g. the transverse DOFs of a pure-axial truss). ``solve_static`` detects them via ``|K_ii| <= 1e-12 * max(|K_ii|)`` and folds them into the Dirichlet set — they drop from the factorisation and return ``u_i = 0`` automatically. This matches MAPDL's default behaviour. See also -------- - :doc:`choosing-a-solver` — which backend for what problem size. - :doc:`performance` — wall-time numbers across backends and meshes. - :doc:`threads` — BLAS / OpenMP thread capping. - :doc:`/user-guide/pre-processing/boundary-conditions` — how to set ``D`` and ``F`` on a Model.