Assembly

FerriteShells.mass_matrix!Function
mass_matrix!(me, scv::ShellCellValues, ρ::T, mat)

Mass matrix for embedded shell elements (2D mesh in 3D). Only translational DOFs contribute to kinetic energy.

source

Kirchhof-Love functions

FerriteShells.membrane_tangent_KL!Function
membrane_tangent_KL!(ke, scv, u_e, mat)

Kirchhoff–Love membrane tangent. u_e is a flat vector of length 3·n_nodes: [$u_1$,$u_2$,$u_3$, $\cdots$].

source
FerriteShells.bending_energy_KLFunction
bending_energy_KL(u_flat, scv::ShellCellValues, mat)

Computes the Kirchhoff–Love bending energy through the curvature change $\kappa_{\alpha\beta} = b_{\alpha\beta} - B_{\alpha\beta}$ (current minus reference second fundamental form).

Arguments:

  • u_flat: a flat vector of length 3·n_nodes containing the displacement degrees of freedom of that element [$u_1$,$u_2$,$u_3$, $\cdots$].
  • scv: a ShellCellValues object containing the precomputed shape function values and derivatives.
  • mat: an instance of an AbstractMaterial describing the material properties.

Note: To capture the full curvature change $\kappa$, quadratic elements are required (linear element only captures twist $\kappa_{12}$).

source

Reissner-Mindlin functions

FerriteShells.membrane_residuals_RM!Function
membrane_residuals_RM!(re, scv, u_e, mat)

RM membrane residual: $r_I = \int N^{\alpha\beta} \partial N_I^\alpha a_\beta \, d\Omega$. Stress resultant rows $P_\alpha = N^{\alpha\beta} a_\beta$ are precomputed once per QP.

source
FerriteShells.bending_residuals_RM!Function
bending_residuals_RM!(re, scv, u_e, mat)

RM bending + transverse shear residual, explicit index-notation form.

Displacement DOFs: $r_I^u = (\partial_1 N_I P^1 + \partial_2 N_I P^2) \, d\Omega$ where $P^\alpha = M^{\alpha\beta} d_{,\beta} + Q^\alpha d$, $M = D : \kappa$ (bending moment), $Q^\alpha = \kappa_s G t A^{\alpha\beta} \gamma_\beta$.

Rotation DOFs: $r_{I,k}^\varphi = F_I \cdot (\partial d_I/\partial \varphi_k) \, d\Omega$ where $F_I = \partial_1 N_I S^1 + \partial_2 N_I S^2 + N_I (Q_1 a_1 + Q_2 a_2)$, $S^\alpha = M^{\alpha\beta} a_\beta$.

source
FerriteShells.bending_shear_energy_RMFunction
bending_shear_energy_RM(u_flat, scv::ShellCellValues, mat)

Reissner–Mindlin bending + transverse shear strain energy. DOF layout: 5 DOFs per node — [$u_1$,$u_2$,$u_3$, $\varphi_1$, $\varphi_2$, $\cdots$] (flat vector of length 5·n_nodes).

Director: $d_I = G_3 + \varphi_{1,I} T_1 + \varphi_{2,I} T_2$ where $G_3$ is the reference unit normal and $T_1$, $T_2$ are reference tangents from scv.

Bending strain: $\kappa_{\alpha\beta} = \frac{1}{2}(a_\alpha \cdot d_{,\beta} + a_\beta \cdot d_{,\alpha}) - B_{\alpha\beta}$

Transverse shear: $\gamma_\alpha = a_\alpha \cdot d$

Shear correction factor $\kappa_s = 5/6$ is applied.

source
FerriteShells.membrane_energy_RMFunction
membrane_energy_RM(u_flat, scv::ShellCellValues, mat)

Reissner–Mindlin membrane strain energy. DOF layout: 5 DOFs per node — [$u_1$,$u_2$,$u_3$, $\varphi_1$, $\varphi_2$, $\cdots$] (flat vector of length 5·n_nodes). Only the displacement DOFs (indices 5I-4:5I-2`) contribute to membrane energy.

source
FerriteShells.membrane_tangent_RM!Function
membrane_tangent_RM!(ke, scv, u_e, mat)

RM membrane tangent. Material part: $K^\text{mat}_{IJ} = \partial N_I^\alpha \partial N_J^\delta M_{\alpha\delta}$ where $M_{\alpha\delta} = C^{\alpha\beta\gamma\delta} a_\beta\otimes a_\gamma$. Geometric part: $K^\text{geo}_{IJ} = (\partial N_I^\alpha N^{\alpha\beta} \partial N_J^\beta) \mathbb{h}_3$. Both $M_{\alpha\delta}$ and $N$ are precomputed once per QP outside the node loops.

source
FerriteShells.bending_tangent_RM!Function
bending_tangent_RM!(ke, scv, u_e, mat)  [MITC dispatch]

Consistent RM bending + transverse shear tangent for MITC elements. The MITC shear sensitivities $\partial\gamma_\alpha/\partial u$ are obtained by differentiating through the tying interpolation:

\[\frac{\partial\gamma_\alpha(q)}{\partial u_J} = \sum_k h^\alpha_{qk} \frac{\partial N_J}{\partial\xi_\alpha}(\xi_k)\cdot d(\xi_k)\]

\[\frac{\partial\gamma_\alpha(q)}{\partial\varphi_{J,l}} = \sum_k h^\alpha_{qk}\, N_J(\xi_k)\left(a_\alpha(\xi_k)\cdot\frac{\partial d_J}{\partial\varphi_l}(\xi_k)\right)\]

where $\partial d_J/\partial\varphi_l(\xi_k)$ is the Rodrigues Jacobian at node J evaluated with the reference geometry $(G_3, T_1, T_2)$ at tying point k. This ensures exact consistency with bending_residuals_RM! so Newton converges quadratically.

Bending ($\kappa$) terms are unchanged from the NoMITC path — only shear ($Q$) terms differ.

source
bending_tangent_RM!(ke, scv, u_e, mat)

RM bending and transverse shear tangent, explicit index-notation form. Four blocks per (I,J) pair:

  • uu (3×3): $\partial_\alpha N_I \partial_\gamma N_J (D^{\alpha\beta\gamma\delta} d_{,\beta} \otimes d_{,\delta}) + q_{IJ}(d \otimes d)$ — frame stiffness with $d_1, d_2$.
  • (3×2): $\partial_\alpha N_I [\delta M^{\alpha\beta} d_{,\beta} + \delta Q^\alpha N_J d] + (g_{IJ}+q_I N_J) \mathrm{dd}_{Jl}$.
  • φu (2×3): filled by transposing the block for (I,J) into the (J,I) position.
  • φφ (2×2): material part $\partial F_I/\partial\varphi_{l,J} \cdot \mathrm{dd}_{Ik}$ + geometric part $F_I \cdot \partial^2 d_I/\partial\varphi_k\partial\varphi_l$ (J=I only).
source

External loading functions

FerriteShells.assemble_traction!Function

Assemble external traction into force vector f for embedded shell elements (2D mesh in 3D). traction is either a Vec{3} (uniform) or a callable x::Vec{3} -> Vec{3}. Uses a FacetQuadratureRule and computes the edge length element directly from 3D node positions, bypassing the sdim mismatch that prevents standard FacetValues from working on embedded meshes. Works for RefQuadrilateral and RefTriangle of any interpolation order.

source
FerriteShells.apply_pointload!Function
apply_pointload!(f, dh, nodeset_name, load)

Add a concentrated force load::Vec{3} to the displacement DOFs of all nodes in nodeset_name. Works for both single-field (:u only) and two-field (:u, ) DofHandlers; in both cases the :u DOFs for node I in a cell occupy local positions 3I-2:3I.

source
FerriteShells.assemble_pressure_tangent!Function
assemble_pressure_tangent!(ke, scv, u_e, p)

Load-stiffness $K_\text{pres} = \partial F_p/\partial u$. Follower pressure $n = a_1 \times a_2$ depends on u through $a_1, a_2$.

\[\frac{\partial n}{\partial u_J} = \partial_1 N_J (e_l \times a_2) + \partial_2 N_J (a_1 \times e_l) = \partial_1 N_J (-\mathrm{spin}(a_2)) + \partial_2 N_J\,\mathrm{spin}(a_1)\]

\[K_{IJ} = p N_I w \bigl[\partial_1 N_J (-\mathrm{spin}(a_2)) + \partial_2 N_J\,\mathrm{spin}(a_1)\bigr]\]

(displacement-displacement block only).

source