sktopt.fea

This package provides functions to handle finite element analysis

class sktopt.fea.FEM_SimpLinearElasticity(task: LinearElasticity, E_min_coeff: float, density_interpolation: Callable = CPUDispatcher(<function simp_interpolation>), solver_option: Literal['spsolve', 'cg_pyamg'] = 'spsolve')

Bases: object

Finite Element solver for linear elasticity using SIMP interpolation.

This class performs linear elastic FEM analysis where the Young’s modulus is interpolated based on material density (ρ) using a SIMP-type interpolation function. It is intended for density-based topology optimization workflows, where element stiffness is expressed as:

E(ρ) = E_min + (E_max - E_min) * f(ρ)

where f(ρ) is typically ρᵖ for SIMP.

Parameters:
  • task (LinearElasticity) – Predefined linear elasticity problem that includes mesh, material constants (E, ν), boundary conditions, load vectors, and basis definitions.

  • E_min_coeff (float) –

    Ratio defining the minimum Young’s modulus as:

    E_min = task.E * E_min_coeff

    Used to avoid singular stiffness matrices during optimization.

  • density_interpolation (Callable, optional) – A function f(ρ) that returns an interpolated stiffness multiplier. Defaults to composer.simp_interpolation (ρᵖ). Any custom interpolation function following SIMP/RAMP/etc. can be used.

  • solver_option ({"spsolve", "cg_pyamg"}, optional) –

    Linear solver backend. - “spsolve”: direct SciPy sparse solver (robust, slower for large DOF) - “cg_pyamg”: Conjugate Gradient with PyAMG multigrid preconditioner

    (fast for large problems)

task

The underlying elasticity problem definition.

Type:

LinearElasticity

E_max

Maximum Young’s modulus (equal to task.E).

Type:

float

E_min

Minimum Young’s modulus used for void regions.

Type:

float

density_interpolation

The SIMP / RAMP interpolation function used to compute material stiffness.

Type:

Callable

solver_option

Selected linear solver backend.

Type:

str

Notes

  • This class does not update densities; it only evaluates the FEM

    response for a given density field.

  • Designed to integrate with OC/MMA/ADMM-based topology optimization

    frameworks.

  • The stiffness matrix assembly depends on interpolated Young’s modulus

    at each element.

  • E_min_coeff should typically be small (1e−3 ~ 1e−9), but not zero.

Examples

>>> fem = FEM_SimpLinearElasticity(
...     task=my_task,
...     E_min_coeff=1e-3,
...     density_interpolation=composer.simp_interpolation,
...     solver_option="spsolve",
... )
>>> u = fem.objectives_multi_load(rho)  # FEM compliaance given density
class sktopt.fea.FEM_SimpLinearHeatConduction(task: LinearHeatConduction, E_min_coeff: float, density_interpolation: Callable = CPUDispatcher(<function simp_interpolation>), solver_option: Literal['spsolve'] = 'spsolve', q: int = 4)

Bases: object

Finite element solver for SIMP-based linear heat conduction problems with support for multi-load objectives and virtual Robin boundaries.

This class evaluates thermal objectives (e.g., average temperature, thermal compliance, or user-defined metrics) for a given density field using SIMP interpolation of the conductivity. It also assembles additional “virtual” Robin boundaries, which enables topology optimization of boundary-dependent heat-transfer behavior.

Parameters:
  • task (LinearHeatConduction) – Problem configuration containing mesh, basis, boundary conditions, load cases, and objective type.

  • E_min_coeff (float) – Minimum conductivity ratio used in the SIMP model. The actual minimum conductivity is task.k * E_min_coeff.

  • density_interpolation (Callable, optional) – Interpolation function for SIMP (or RAMP) mapping ρ → k(ρ). Defaults to composer.simp_interpolation.

  • solver_option ({"spsolve", "cg_pyamg"}, optional) – Linear solver to use for the state equation of each load case.

  • q (int, optional) – Exponent for boundary interpolation in the virtual Robin model (often used to sharpen on/off behavior of boundary heat transfer).

k_max

Conductivity of solid material.

Type:

float

k_min

Conductivity of void (or weak) material.

Type:

float

λ_all

Stored adjoint fields for all load cases, computed during the last call to objectives_multi_load().

Type:

np.ndarray or None

Notes

  • The class converts element-wise densities into nodal densities by averaging over connected elements. This ensures smoother interpolation when assembling virtual Robin terms.

  • When task.design_robin_boundary is True, the actual boundary condition is updated based on the given density field, enabling optimization over Robin boundaries.

  • For multi-load problems, each load case is solved independently and the objectives are returned as a list (e.g., one compliance value per thermal load case).