QUAD4_PLANE — 2-D 4-node structural solid#
Kinematics. Bilinear 2-D quadrilateral; two translational DOFs per node (\(u_x, u_y\)); 8 DOFs per element. Lives in the global (x, y) plane.
Stiffness. \(K_e = t \int_A B^\top C\, B\,\mathrm{d}A\) integrated on the reference square with 2 × 2 Gauss.
Mass. Consistent \(M_e = \rho\, t \int_A N^\top N\,\mathrm{d}A\) with 2 × 2 Gauss.
Modes. Plane stress (default) and plane strain are implemented today; the axisymmetric mode is roadmap.
Bilinear shape functions#
On the reference square \((\xi, \eta) \in [-1, 1]^2\), the four corner nodes at \((\xi_i, \eta_i) \in \{\pm 1\}^2\) in CCW order carry the bilinear shape functions [Cook2001] [ZTZ2013]
Node order matches VTK_QUAD — the same \((-,-), (+,-), (+,+), (-,+)\) sweep used by QUAD4_SHELL.
Jacobian, \(B\) matrix, and Voigt strain#
Two-dimensional small strain with engineering shear: \(\varepsilon = [\varepsilon_{xx},\ \varepsilon_{yy},\ \gamma_{xy}]^\top\). The 3 × 8 strain-displacement matrix assembles from two-row blocks per node:
Chain-rule derivatives come from the 2 × 2 Jacobian \(J_{jk} = \partial x_k / \partial \xi_j\).
Plane-stress vs plane-strain constitutive#
The mode argument on ELEMENTS.QUAD4_PLANE selects:
Plane stress (
mode="stress", default): thin in \(\hat{z}\), \(\sigma_{zz} = 0\).\[\begin{split}C_\text{ps} = \frac{E}{1 - \nu^2} \begin{bmatrix} 1 & \nu & 0 \\ \nu & 1 & 0 \\ 0 & 0 & (1 - \nu) / 2 \end{bmatrix}.\end{split}\]Plane strain (
mode="strain"): long in \(\hat{z}\), \(\varepsilon_{zz} = 0\).\[\begin{split}C_\text{pe} = \frac{E (1 - \nu)}{(1 + \nu)(1 - 2\nu)} \begin{bmatrix} 1 & \tfrac{\nu}{1 - \nu} & 0 \\ \tfrac{\nu}{1 - \nu} & 1 & 0 \\ 0 & 0 & \tfrac{1 - 2\nu}{2 (1 - \nu)} \end{bmatrix}.\end{split}\]
2 × 2 Gauss quadrature#
The 2 × 2 rule with points at \((\pm 1/\sqrt{3}, \pm 1/\sqrt{3})\) and unit weights integrates bilinear \(B^\top C B\) exactly on undistorted quads.
Thickness \(t\) comes from the element’s real-constant slot
real[0] (default 1.0 for plane-strain use).
Numpy walk-through#
import numpy as np
xi_n = np.array([(-1,-1), (+1,-1), (+1,+1), (-1,+1)], dtype=float)
def shape_and_derivs(xi, eta):
N = 0.25 * (1 + xi_n[:, 0] * xi) * (1 + xi_n[:, 1] * eta)
dN = np.empty((4, 2))
dN[:, 0] = 0.25 * xi_n[:, 0] * (1 + xi_n[:, 1] * eta)
dN[:, 1] = 0.25 * xi_n[:, 1] * (1 + xi_n[:, 0] * xi)
return N, dN
g = 1.0 / np.sqrt(3.0)
gp = np.array([(-g,-g), (+g,-g), (+g,+g), (-g,+g)])
w = np.ones(4)
# Plane stress
E, nu = 2.1e11, 0.30
C = (E / (1 - nu**2)) * np.array([[1, nu, 0], [nu, 1, 0], [0, 0, (1 - nu) / 2]])
# Unit-square nodes in the (x, y) plane.
X = np.array([(0, 0), (1, 0), (1, 1), (0, 1)], dtype=float)
K = np.zeros((8, 8))
thick = 1.0
for (xi, eta), wi in zip(gp, w):
_, dNdxi = shape_and_derivs(xi, eta)
J = dNdxi.T @ X
detJ = np.linalg.det(J)
dNdx = np.linalg.solve(J, dNdxi.T).T
B = np.zeros((3, 8))
for i in range(4):
B[0, 2 * i ] = dNdx[i, 0]
B[1, 2 * i + 1] = dNdx[i, 1]
B[2, 2 * i ] = dNdx[i, 1]
B[2, 2 * i + 1] = dNdx[i, 0]
K += thick * (B.T @ C @ B) * detJ * wi
\(K_e\) has 3 rigid-body zeros (two translations + one rotation in the plane) and 5 elastic modes on an undistorted unit square.
Mode flag#
|
Description |
|---|---|
|
Plane stress; thin in \(\hat{z}\). Default thickness = 1.0. |
|
Plane strain; long in \(\hat{z}\), \(\varepsilon_{zz} = 0\). |
The axisymmetric mode is not yet in the present kernel.
Validation#
Rigid-body modes. 3 rigid-body zero eigenvalues + 5 elastic modes on an undistorted square.
Patch test. A linear displacement field ramps exactly through
every node of any distorted quad mesh — verified in
tests.elements.plane182.
Convergence. Mesh refinement drives in-plane cantilever
deflection to the plane-stress closed form under both
"plane_stress" and "plane_strain" (with the appropriate
analytic reference per mode).
API reference#
- class femorph_solver.elements.quad4_plane.Quad4Plane[source]#
Bases:
ElementBase- static eel_batch(coords: ndarray, u_e: ndarray, material: dict[str, float] | None = None) ndarray | None[source]#
Per-element engineering-shear Voigt strain at element nodes.
coords:(n_elem, 4, 3)or(n_elem, 4, 2).u_e:(n_elem, 8)— DOF layout[ux0, uy0, ux1, …].Returns
(n_elem, 4, 6)Voigt strain ordered[εxx, εyy, εzz, γxy, γyz, γxz]— matches the 3D-solid Voigt convention socompute_nodal_stress()can drive the sameσ = C·εcontraction with the standard isotropic(6, 6)C.For plane stress (default) we fill
εzz = -ν/(1-ν)·(εxx+εyy)soσ_zzrecovered via the 6×6Cis identically zero (the defining condition of plane stress). For plane strain we fillεzz = 0(the defining condition of plane strain);σ_zzis thenν·(σ_xx + σ_yy)per Hooke’s law.materialmust be provided when the kernel is operating in plane stress (the default) soνis available; for plane strain thematerialargument is unused.Strain is evaluated at each corner node (ξ_i, η_i) via the bilinear B-matrix — the same basis used in
ke. This is the standard Q4 nodal recovery; the recovered strain is exact at the centroid and accurate toO(h²)at corner nodes for smooth solutions.
References#
Cook, R. D., Malkus, D. S., Plesha, M. E., and Witt, R. J. (2001). Concepts and Applications of Finite Element Analysis, 4th ed. Wiley. Ch. 3 (plane stress / plane strain). https://www.wiley.com/en-us/Concepts+and+Applications+of+Finite+Element+Analysis%2C+4th+Edition-p-9780471356059
Zienkiewicz, O. C., Taylor, R. L., and Zhu, J. Z. (2013). The Finite Element Method: Its Basis and Fundamentals, 7th ed. Butterworth-Heinemann. Ch. 5 (isoparametric plane elements). https://doi.org/10.1016/C2009-0-24909-9
Bathe, K.-J. (2014). Finite Element Procedures, 2nd ed. https://www.klausjurgenbathe.com/fepbook/