fenics_shells.common package¶
Submodules¶
fenics_shells.common.constitutive_models module¶
-
fenics_shells.common.constitutive_models.
psi_M
(k, **kwargs)¶ Returns bending moment energy density calculated from the curvature k using:
Isotropic case: .. math:: D = frac{E*t^3}{24(1 - nu^2)} W_m(k, ldots) = D*((1 - nu)*tr(k**2) + nu*(tr(k))**2)
Parameters: - k – Curvature, typically UFL form with shape (2,2) (tensor).
- **kwargs – Isotropic case: E: Young’s modulus, Constant or Expression. nu: Poisson’s ratio, Constant or Expression. t: Thickness, Constant or Expression.
Returns: UFL form of bending stress tensor with shape (2,2) (tensor).
-
fenics_shells.common.constitutive_models.
psi_N
(e, **kwargs)¶ Returns membrane energy density calculated from e using:
Isotropic case: .. math:: B = frac{E*t}{2(1 - nu^2)} N(e, ldots) = B(1 - nu)e + nu mathrm{tr}(e)I
Parameters: - e – Membrane strain, typically UFL form with shape (2,2) (tensor).
- **kwargs – Isotropic case: E: Young’s modulus, Constant or Expression. nu: Poisson’s ratio, Constant or Expression. t: Thickness, Constant or Expression.
Returns: UFL form of membrane stress tensor with shape (2,2) (tensor).
-
fenics_shells.common.constitutive_models.
strain_from_voigt
(e_voigt)¶ Inverse operation of strain_to_voigt.
Parameters: - sigma_voigt – UFL form with shape (3,1) corresponding to the strain
- in Voigt format (pseudo-vector) –
Returns: a symmetric stress tensor, typically UFL form with shape (2,2)
-
fenics_shells.common.constitutive_models.
strain_to_voigt
(e)¶ Returns the pseudo-vector in the Voigt notation associate to a 2x2 symmetric strain tensor, according to the following rule (see e.g. https://en.wikipedia.org/wiki/Voigt_notation),
\[\begin{split}e = \begin{bmatrix} e_{00} & e_{01}\\ e_{01} & e_{11} \end{bmatrix}\quad\to\quad e_\mathrm{voigt}= \begin{bmatrix} e_{00} & e_{11}& 2e_{01} \end{bmatrix}\end{split}\]Parameters: e – a symmetric 2x2 strain tensor, typically UFL form with shape (2,2) Returns: a UFL form with shape (3,1) corresponding to the input tensor in Voigt notation.
-
fenics_shells.common.constitutive_models.
stress_from_voigt
(sigma_voigt)¶ Inverse operation of stress_to_voigt.
Parameters: - sigma_voigt – UFL form with shape (3,1) corresponding to the stress
- in Voigt format. (pseudo-vector) –
Returns: a symmetric stress tensor, typically UFL form with shape (2,2)
-
fenics_shells.common.constitutive_models.
stress_to_voigt
(sigma)¶ Returns the pseudo-vector in the Voigt notation associate to a 2x2 symmetric stress tensor, according to the following rule (see e.g. https://en.wikipedia.org/wiki/Voigt_notation),
\[\begin{split}\sigma = \begin{bmatrix} \sigma_{00} & \sigma_{01}\\ \sigma_{01} & \sigma_{11} \end{bmatrix}\quad\to\quad \sigma_\mathrm{voigt}= \begin{bmatrix} \sigma_{00} & \sigma_{11}& \sigma_{01} \end{bmatrix}\end{split}\]Parameters: - sigma – a symmetric 2x2 stress tensor, typically UFL form with shape
- (2,2) –
Returns: a UFL form with shape (3,1) corresponding to the input tensor in Voigt notation.
fenics_shells.common.energy module¶
-
fenics_shells.common.energy.
membrane_bending_energy
(e, k, A, D, B)¶ Return the coupled membrane-bending energy for a plate model.
Parameters: - e – Membrane strains, UFL or DOLFIN Function of rank (2, 2) (tensor).
- k – Curvature, UFL or DOLFIN Function of rank (2, 2) (tensor).
- A – Membrane stresses.
- D – Bending stresses.
- B – Coupled membrane-bending stresses.
-
fenics_shells.common.energy.
membrane_energy
(e, N)¶ Return internal membrane energy for a plate model.
Parameters: - e – Membrane strains, UFL or DOLFIN Function of rank (2, 2) (tensor).
- N – Membrane stress, UFL or DOLFIN Function of rank (2, 2) (tensor).
Returns: UFL form of internal elastic membrane energy for a plate model.
fenics_shells.common.kinematics module¶
-
fenics_shells.common.kinematics.
F
(u)¶ Return deformation gradient tensor for non-linear plate model.
Deformation gradient of 2-dimensional manifold embedded in 3-dimensional space.
\[F = I + \nabla u\]Parameters: u – displacement field, typically UFL (3,1) coefficient. Returns: a UFL coeffient with shape (3,2)
-
fenics_shells.common.kinematics.
e
(u)¶ Return membrane strain tensor for linear plate model.
\[e = \dfrac{1}{2}(\nabla u+\nabla u^T)\]Parameters: u – membrane displacement field, typically UFL (2,1) coefficient. Returns: a UFL form with shape (2,2)
-
fenics_shells.common.kinematics.
k
(theta)¶ Return bending curvature tensor for linear plate model.
\[k = \dfrac{1}{2}(\nabla \theta+\nabla \theta^T)\]Parameters: theta – rotation field, typically UFL (2,1) form or a dolfin Function Returns: a UFL form with shape (2,2)
fenics_shells.common.laminates module¶
-
fenics_shells.common.laminates.
ABD
(E1, E2, G12, nu12, hs, thetas)¶ Return the stiffness matrix of a kirchhoff-love model of a laminate obtained by stacking n orthotropic laminae with possibly different thinknesses and orientations (see Reddy 1997, eqn 1.3.71).
It assumes a plane-stress state.
Parameters: - E1 – The Young modulus in the material direction 1.
- E2 – The Young modulus in the material direction 2.
- G12 – The in-plane shear modulus.
- nu12 – The in-plane Poisson ratio.
- hs – a list with length n with the thicknesses of the layers (from top to bottom).
- theta – a list with the n orientations (in radians) of the layers (from top to bottom).
Returns: a symmetric 3x3 ufl matrix giving the membrane stiffness in Voigt notation. B: a symmetric 3x3 ufl matrix giving the membrane/bending coupling stiffness in Voigt notation. D: a symmetric 3x3 ufl matrix giving the bending stiffness in Voigt notation.
Return type: A
-
fenics_shells.common.laminates.
F
(G13, G23, hs, thetas)¶ Return the shear stiffness matrix of a Reissner-Midlin model of a laminate obtained by stacking n orthotropic laminae with possibly different thinknesses and orientations. (See Reddy 1997, eqn 3.4.18)
It assumes a plane-stress state.
Parameters: - G13 – The transverse shear modulus between the material directions 1-3.
- G23 – The transverse shear modulus between the material directions 2-3.
- hs – a list with length n with the thicknesses of the layers (from top to bottom).
- theta – a list with the n orientations (in radians) of the layers (from top to bottom).
Returns: a symmetric 2x2 ufl matrix giving the shear stiffness in Voigt notation.
Return type:
-
fenics_shells.common.laminates.
NM_T
(E1, E2, G12, nu12, hs, thetas, DeltaT_0, DeltaT_1=0.0, alpha1=1.0, alpha2=1.0)¶ Return the thermal stress and moment resultant of a Kirchhoff-Love model of a laminate obtained by stacking n orthotropic laminae with possibly different thinknesses and orientations.
It assumes a plane-stress states and a temperature distribution in the from
Delta(z) = DeltaT_0 + z * DeltaT_1
Parameters: - E1 – The Young modulus in the material direction 1.
- E2 – The Young modulus in the material direction 2.
- G12 – The in-plane shear modulus.
- nu12 – The in-plane Poisson ratio.
- hs – a list with length n with the thicknesses of the layers (from top to bottom).
- theta – a list with the n orientations (in radians) of the layers (from top to bottom).
- alpha1 – Expansion coefficient in the material direction 1.
- alpha2 – Expansion coefficient in the material direction 2.
- DeltaT_0 – Average temperature field.
- DeltaT_1 – Gradient of the temperature field.
Returns: a 3x1 ufl vector giving the membrane inelastic stress. M_T: a 3x1 ufl vector giving the bending inelastic stress.
Return type: N_T
-
fenics_shells.common.laminates.
rotated_lamina_expansion_inplane
(alpha11, alpha22, theta)¶ Return the in-plane expansion matrix of an orhtropic layer in a reference rotated by an angle theta wrt to the material one. It assumes Voigt notation and plane stress state. (See Reddy 1997, eqn 1.3.71)
Parameters: - alpha11 – Expansion coefficient in the material direction 1.
- alpha22 – Expansion coefficient in the material direction 2.
- theta – The rotation angle from the material to the desired reference system.
Returns: a 3x1 ufl vector giving the expansion matrix in voigt notation.
Return type: alpha_theta
-
fenics_shells.common.laminates.
rotated_lamina_stiffness_inplane
(E1, E2, G12, nu12, theta)¶ Return the in-plane stiffness matrix of an orhtropic layer in a reference rotated by an angle theta wrt to the material one. It assumes Voigt notation and plane stress state. (See Reddy 1997, eqn 1.3.71)
Parameters: - E1 – The Young modulus in the material direction 1.
- E2 – The Young modulus in the material direction 2.
- G23 – The in-plane shear modulus
- nu12 – The in-plane Poisson ratio
- theta – The rotation angle from the material to the desired refence system
Returns: a 3x3 symmetric ufl matrix giving the stiffness matrix
Return type: Q_theta
-
fenics_shells.common.laminates.
rotated_lamina_stiffness_shear
(G13, G23, theta, kappa=0.8333333333333334)¶ Return the shear stiffness matrix of an orhtropic layer in a reference rotated by an angle theta wrt to the material one. It assumes Voigt notation and plane stress state (see Reddy 1997, eqn 3.4.18).
Parameters: - G12 – The transverse shear modulus between the material directions 1-2.
- G13 – The transverse shear modulus between the material directions 1-3.
- kappa – The shear correction factor.
Returns: a 3x3 symmetric ufl matrix giving the stiffness matrix.
Return type: Q_shear_theta
-
fenics_shells.common.laminates.
z_coordinates
(hs)¶ Return a list with the thickness coordinate of the top surface of each layer taking the midplane as z = 0.
Parameters: hs – a list giving the thinckesses of each layer ordered from bottom (layer - 0) to top (layer n-1). Returns: - a list of coordinate of the top surface of each layer
- ordered from bottom (layer - 0) to top (layer n-1)
Return type: z