Physics#
Physics-level helpers (constitutive models, material laws, etc.).
- fluxfem.physics.assemble_constant_body_force(space: FESpaceClosure, gravity_vec: jax.numpy.ndarray, density: float, *, sparse: bool = False) jax.numpy.ndarray | tuple[jax.numpy.ndarray, jax.numpy.ndarray, int]#
Convenience: assemble body force from density * gravity vector. gravity_vec: length-3 array-like (direction and magnitude of g) density: scalar density (consistent with unit system)
- fluxfem.physics.constant_body_force_vector_form(ctx: jax.tree_util.register_pytree_node_class, load_vec: jax.numpy.ndarray) jax.numpy.ndarray#
Linear form for 3D vector body force f (constant in space).
- fluxfem.physics.ddot(a: ArrayLike, b: ArrayLike, c: ArrayLike | None = None) ArrayLike#
Double contraction on the last two axes.
ddot(a, b): sum_ij a_ij * b_ij
ddot(a, b, c): a^T b c (Voigt-style linear elasticity blocks)
- fluxfem.physics.deformation_gradient(ctx: jax.tree_util.register_pytree_node_class, u_elem: jax.numpy.ndarray) jax.numpy.ndarray#
Compute deformation gradient F = I + grad_u for a 3D vector displacement.
ctx: FormContext with test/trial grads (reference configuration) u_elem: (n_ldofs,) element displacement in dof ordering [u0,v0,w0, u1,v1,w1, …] returns: (n_q, 3, 3) F per quadrature point
- fluxfem.physics.diffusion_form(ctx: jax.tree_util.register_pytree_node_class, kappa: float) jax.numpy.ndarray#
Scalar diffusion bilinear form: kappa * grad_v · grad_u.
Returns the per-quadrature integrand for a standard Laplacian term.
- fluxfem.physics.dot(a: FormFieldLike | ArrayLike, b: ArrayLike) ArrayLike#
Batched matrix product on the last two axes.
If the first argument is a FormField, dispatch to vector_load_form to build the linear form contribution for a vector load.
- fluxfem.physics.green_lagrange_strain(F: jax.numpy.ndarray) jax.numpy.ndarray#
E = 0.5 (C - I).
- fluxfem.physics.interpolate_at_points(space: FESpaceClosure, u: ndarray, points: ndarray) ndarray#
Interpolate displacement field at given physical points (Hex8 only, structured search). - points: (m,3) array of physical coordinates. Returns: (m,3) interpolated displacement.
- fluxfem.physics.isotropic_3d_D(E: float, nu: float) jax.numpy.ndarray#
Return 3D isotropic linear elasticity constitutive matrix in Voigt form.
- fluxfem.physics.lame_parameters(E: float, nu: float) tuple[float, float]#
Return Lamé parameters (lambda, mu) from Young’s modulus and Poisson ratio.
- fluxfem.physics.linear_elasticity_form(ctx: jax.tree_util.register_pytree_node_class, D: jax.numpy.ndarray) jax.numpy.ndarray#
Linear-elasticity bilinear form in Voigt notation.
Returns the per-quadrature integrand for Bv^T D Bu, where B is the symmetric-gradient operator for the test/trial fields.
- fluxfem.physics.make_elastic_point_data(mesh: BaseMeshClosure, space: FESpaceClosure, u: ndarray, *, compute_j: bool = True, deformed_scale: float = 1.0) dict[str, ndarray]#
Alias to postprocess.make_point_data_displacement for backward compatibility.
- fluxfem.physics.make_point_data_displacement(mesh: BaseMeshClosure, space: FESpaceClosure, u: ndarray, *, compute_j: bool = True, deformed_scale: float = 1.0) dict[str, ndarray]#
- Common postprocess helper to build point data dictionaries:
displacement
deformed_coords = X + scale * u
optional J (nodal average of det(F))
Assumes 3 dof/node ordering [u0,v0,w0, u1,v1,w1, …].
- fluxfem.physics.max_shear_stress(S: jax.numpy.ndarray) jax.numpy.ndarray#
Maximum shear stress = (sigma_max - sigma_min) / 2.
- fluxfem.physics.neo_hookean_residual_form(ctx: jax.tree_util.register_pytree_node_class, u_elem: jax.numpy.ndarray, params: Mapping[str, float] | tuple[float, float]) jax.numpy.ndarray#
Compressible Neo-Hookean residual (Total Lagrangian, PK2). params: dict-like with keys “mu”, “lam” or tuple (mu, lam) returns: (n_q, n_ldofs)
- fluxfem.physics.pk2_neo_hookean(F: jax.numpy.ndarray, mu: float, lam: float) jax.numpy.ndarray#
- Compressible Neo-Hookean PK2 stress:
S = mu * (I - C^{-1}) + lam * ln J * C^{-1} C = F^T F, J = sqrt(det C)
- fluxfem.physics.principal_stresses(S: jax.numpy.ndarray) jax.numpy.ndarray#
Return principal stresses (eigvals) for symmetric 3x3 stress tensor. Supports batching over leading dimensions.
- fluxfem.physics.principal_sum(S: jax.numpy.ndarray) jax.numpy.ndarray#
Sum of principal stresses (trace).
- fluxfem.physics.right_cauchy_green(F: jax.numpy.ndarray) jax.numpy.ndarray#
C = F^T F (right Cauchy-Green).
- fluxfem.physics.sym_grad(field: FormFieldLike) jax.numpy.ndarray#
Symmetric gradient operator for vector mechanics (small strain).
- Parameters:
field (FormField-like) –
- Must provide:
field.gradN : (n_q, n_nodes, 3)
field.basis.dofs_per_node (usually 3)
- Returns:
B – (n_q, 6, dofs_per_node*n_nodes) Voigt order [xx, yy, zz, xy, yz, zx] Such that eps_voigt(q,:) = B(q,:,:) @ u_elem
- Return type:
jnp.ndarray
- fluxfem.physics.sym_grad_u(field: FormFieldLike, u_elem: jax.numpy.ndarray) jax.numpy.ndarray#
Apply sym_grad(field) to a local displacement vector.
- Parameters:
field (FormField-like) – Vector field basis data.
u_elem (jnp.ndarray) – Element displacement vector (dofs_per_node*n_nodes,).
- Returns:
Symmetric strain in Voigt form with shape (n_q, 6).
- Return type:
jnp.ndarray
- fluxfem.physics.transpose_last2(a: jax.numpy.ndarray) jax.numpy.ndarray#
Swap the last two axes (batched transpose).
- fluxfem.physics.vector_body_force_form(ctx: jax.tree_util.register_pytree_node_class, load_vec: jax.numpy.ndarray) jax.numpy.ndarray#
Linear form for 3D vector body force f (constant in space).
- fluxfem.physics.von_mises_stress(S: jax.numpy.ndarray) jax.numpy.ndarray#
von Mises equivalent stress: sqrt(3/2 * dev(S):dev(S)).
- fluxfem.physics.write_elastic_vtu(mesh: BaseMeshClosure, space: FESpaceClosure, u: ndarray, filepath: str, *, compute_j: bool = True, deformed_scale: float = 1.0) None#
Alias to postprocess.write_point_data_vtu for backward compatibility.
- fluxfem.physics.write_point_data_vtu(mesh: BaseMeshClosure, space: FESpaceClosure, u: ndarray, filepath: str, *, compute_j: bool = True, deformed_scale: float = 1.0) None#
Write VTU with displacement/deformed_coords and optional J.
Diffusion#
- fluxfem.physics.diffusion.diffusion_form(ctx: jax.tree_util.register_pytree_node_class, kappa: float) jax.numpy.ndarray#
Scalar diffusion bilinear form: kappa * grad_v · grad_u.
Returns the per-quadrature integrand for a standard Laplacian term.
Operators#
- fluxfem.physics.operators.ddot(a: ArrayLike, b: ArrayLike, c: ArrayLike | None = None) ArrayLike#
Double contraction on the last two axes.
ddot(a, b): sum_ij a_ij * b_ij
ddot(a, b, c): a^T b c (Voigt-style linear elasticity blocks)
- fluxfem.physics.operators.dot(a: FormFieldLike | ArrayLike, b: ArrayLike) ArrayLike#
Batched matrix product on the last two axes.
If the first argument is a FormField, dispatch to vector_load_form to build the linear form contribution for a vector load.
- fluxfem.physics.operators.sym_grad(field: FormFieldLike) jax.numpy.ndarray#
Symmetric gradient operator for vector mechanics (small strain).
- Parameters:
field (FormField-like) –
- Must provide:
field.gradN : (n_q, n_nodes, 3)
field.basis.dofs_per_node (usually 3)
- Returns:
B – (n_q, 6, dofs_per_node*n_nodes) Voigt order [xx, yy, zz, xy, yz, zx] Such that eps_voigt(q,:) = B(q,:,:) @ u_elem
- Return type:
jnp.ndarray
- fluxfem.physics.operators.sym_grad_u(field: FormFieldLike, u_elem: jax.numpy.ndarray) jax.numpy.ndarray#
Apply sym_grad(field) to a local displacement vector.
- Parameters:
field (FormField-like) – Vector field basis data.
u_elem (jnp.ndarray) – Element displacement vector (dofs_per_node*n_nodes,).
- Returns:
Symmetric strain in Voigt form with shape (n_q, 6).
- Return type:
jnp.ndarray
- fluxfem.physics.operators.transpose_last2(a: jax.numpy.ndarray) jax.numpy.ndarray#
Swap the last two axes (batched transpose).
Postprocess#
- fluxfem.physics.postprocess.interpolate_at_points(space: FESpaceClosure, u: ndarray, points: ndarray) ndarray#
Interpolate displacement field at given physical points (Hex8 only, structured search). - points: (m,3) array of physical coordinates. Returns: (m,3) interpolated displacement.
- fluxfem.physics.postprocess.make_point_data_displacement(mesh: BaseMeshClosure, space: FESpaceClosure, u: ndarray, *, compute_j: bool = True, deformed_scale: float = 1.0) dict[str, ndarray]#
- Common postprocess helper to build point data dictionaries:
displacement
deformed_coords = X + scale * u
optional J (nodal average of det(F))
Assumes 3 dof/node ordering [u0,v0,w0, u1,v1,w1, …].
- fluxfem.physics.postprocess.write_point_data_vtu(mesh: BaseMeshClosure, space: FESpaceClosure, u: ndarray, filepath: str, *, compute_j: bool = True, deformed_scale: float = 1.0) None#
Write VTU with displacement/deformed_coords and optional J.