Hyperelastic Plane Stress¶
Bases: Hyperelastic3D
Hyperelastic plane stress material.
This class extends Hyperelastic3D for plane stress problems by enforcing
the out-of-plane stress condition \(P_{33} = 0\) through a local Newton
iteration on the out-of-plane stretch \(\lambda_3\).
Parameters:
-
psi(Callable) –Strain energy density function \(\psi(\mathbf{F}, \texttt{params})\). This must be a 3D energy; the plane stress condensation is performed internally.
-
params(list | Tensor) –Material parameters passed to
psi. Shape:(p,)for a scalar or(N, p)for a batch of materials. -
rho(Tensor | float, default:1.0) –Mass density. Default is
1.0. -
tolerance(float, default:1e-05) –Convergence tolerance for the local Newton solver. Default is
1e-5. -
max_iter(int, default:10) –Maximum number of iterations for the local Newton solver. Default is
10.
Notes
- Finite-strain (large deformation) framework.
- One internal state variable (
n_state = 1): the out-of-plane stretch increment \(\lambda_3 - 1\). - The 3D deformation gradient is reconstructed as \(\mathbf{F} = \text{diag}(\mathbf{F}^{2D},\, \lambda_3)\).
Plane stress condensation
Given a 3D strain energy \(\psi(\mathbf{F})\), the plane stress condition requires \(P_{33} = 0\). The out-of-plane stretch \(\lambda_3\) is found iteratively by solving $$ P_{33} (\mathbf{F}^{\text{2D}}, \lambda_3) = 0 $$ with a Newton update $$ \lambda_3 \leftarrow \lambda_3 - \frac{P_{33}}{C_{3333}}. $$ The in-plane tangent is obtained via static condensation: $$ C_{ijkl}^{\text{2D}} = C_{ijkl} - \frac{C_{ij33}\,C_{33kl}}{C_{3333}}. $$
vectorize(n_elem)
¶
Returns a vectorized copy of the material for n_elem elements.
This function creates a batched version of the material properties. If the
material is already vectorized (self.is_vectorized == True), the function
simply returns self without modification.
Parameters:
-
n_elem(int) –Number of elements to vectorize the material for.
Returns:
-
HyperelasticPlaneStress(HyperelasticPlaneStress) –A new material instance with vectorized properties.
step(H_inc, F, stress, state, de0, cl, iter)
¶
Performs an incremental step for a hyperelastic material.
Parameters:
-
H_inc(Tensor) –Incremental displacement gradient. - Shape:
(..., 2, 2), where...represents batch dimensions. -
F(Tensor) –Current deformation gradient. - Shape:
(..., 2, 2), same asH_inc. -
stress(Tensor) –Current 1st PK stress tensor. - Shape:
(..., 2, 2). -
state(Tensor) –Internal state variables (unused in hyperelasticity). - Shape: Arbitrary, remains unchanged.
-
de0(Tensor) –External deformation gradient increment (e.g., thermal). - Shape:
(..., 2, 2). -
cl(Tensor) –Characteristic lengths. - Shape:
(..., 1). -
iter(int) –Current iteration number.
Returns:
-
tuple(tuple[Tensor, Tensor, Tensor]) –- stress_new (Tensor): Updated 1st PK stress tensor.
Shape:
(..., 2, 2). - state_new (Tensor): Updated internal state (unchanged).
Shape: same as
state. - ddsdde (Tensor): Algorithmic tangent stiffness tensor.
Shape:
(..., 2, 2, 2, 2).
- stress_new (Tensor): Updated 1st PK stress tensor.
Shape: