sfepy.mechanics.membranes module¶
-
sfepy.mechanics.membranes.create_mapping(coors, gel, order)[source]¶ Create mapping from transformed (in x-y plane) element faces to reference element faces.
Parameters: coors : array
The transformed coordinates of element nodes, shape (n_el, n_ep, dim). The function verifies that the all z components are zero.
gel : GeometryElement instance
The geometry element corresponding to the faces.
order : int
The polynomial order of the mapping.
Returns: mapping : VolumeMapping instance
The reference element face mapping.
-
sfepy.mechanics.membranes.create_transformation_matrix(coors)[source]¶ Create a transposed coordinate transformation matrix, that transforms 3D coordinates of element face nodes so that the transformed nodes are in the x-y plane. The rotation is performed w.r.t. the first node of each face.
Parameters: coors : array
The coordinates of element nodes, shape (n_el, n_ep, dim).
Returns: mtx_t : array
The transposed transformation matrix
, i.e.
.Notes
, where
,
, are unit
in-plane (column) vectors and
is the unit normal vector,
all mutually orthonormal.
-
sfepy.mechanics.membranes.describe_deformation(el_disps, bfg)[source]¶ Describe deformation of a thin incompressible 2D membrane in 3D space, composed of flat finite element faces.
The coordinate system of each element (face), i.e. the membrane mid-surface, should coincide with the x, y axes of the x-y plane.
Parameters: el_disps : array
The displacements of element nodes, shape (n_el, n_ep, dim).
bfg : array
The in-plane base function gradients, shape (n_el, n_qp, dim-1, n_ep).
Returns: mtx_c ; array
The in-plane right Cauchy-Green deformation tensor
,
.c33 : array
The component
computed from the incompressibility
condition.mtx_b : array
The discrete Green strain variation operator.
-
sfepy.mechanics.membranes.describe_geometry(ig, field, region, integral)[source]¶ Describe membrane geometry in a given region.
Parameters: ig : int
The element group index.
field : Field instance
The field defining the FE approximation.
region : Region instance
The surface region to describe.
integral : Integral instance
The integral defining the quadrature points.
Returns: mtx_t : array
The transposed transformation matrix
, see
create_transformation_matrix().membrane_geo : CMapping instance
The mapping from transformed elements to a reference elements.
-
sfepy.mechanics.membranes.get_green_strain_sym3d(mtx_c, c33)[source]¶ Get the 3D Green strain tensor in symmetric storage.
Parameters: mtx_c ; array
The in-plane right Cauchy-Green deformation tensor
,
, shape (n_el, n_qp, dim-1,
dim-1).c33 : array
The component
computed from the incompressibility
condition, shape (n_el, n_qp).Returns: mtx_e : array
The membrane Green strain
, symmetric storage: items (11, 22, 33, 12, 13, 23),
shape (n_el, n_qp, sym, 1).
-
sfepy.mechanics.membranes.get_invariants(mtx_c, c33)[source]¶ Get the first and second invariants of the right Cauchy-Green deformation tensor describing deformation of an incompressible membrane.
Parameters: mtx_c ; array
The in-plane right Cauchy-Green deformation tensor
,
, shape (n_el, n_qp, dim-1,
dim-1).c33 : array
The component
computed from the incompressibility
condition, shape (n_el, n_qp).Returns: i1 : array
The first invariant of
.i2 : array
The second invariant of
.
-
sfepy.mechanics.membranes.get_tangent_stress_matrix(stress, bfg)[source]¶ Get the tangent stress matrix of a thin incompressible 2D membrane in 3D space, given a stress.
Parameters: stress : array
The components 11, 22, 12 of the second Piola-Kirchhoff stress tensor, shape (n_el, n_qp, 3, 1).
bfg : array
The in-plane base function gradients, shape (n_el, n_qp, dim-1, n_ep).
Returns: mtx : array
The tangent stress matrix, shape (n_el, n_qp, dim*n_ep, dim*n_ep).

