fenics_shells.reissner_mindlin package

Submodules

fenics_shells.reissner_mindlin.forms module

fenics_shells.reissner_mindlin.forms.gamma(theta, w)

Return shear strain vector calculated from primal variables:

\[\gamma = \nabla w - \theta\]
fenics_shells.reissner_mindlin.forms.inner_e(x, y, restrict_to_one_side=False, quadrature_degree=1)

The inner product of the tangential component of a vector field on all of the facets of the mesh (Measure objects dS and ds).

By default, restrict_to_one_side is False. In this case, the function will return an integral that is restricted to both sides (‘+’) and (‘-’) of a shared facet between elements. You should use this in the case that you want to use the ‘projected’ version of DuranLibermanSpace.

If restrict_to_one_side is True, then this will return an integral that is restricted (‘+’) to one side of a shared facet between elements. You should use this in the case that you want to use the multipliers version of DuranLibermanSpace.

Parameters:
  • x – DOLFIN or UFL Function of rank (2,) (vector).
  • y – DOLFIN or UFL Function of rank (2,) (vector).
  • (Optional[bool] (restrict_to_one_side) – Default is False.
  • quadrature_degree (Optional[int]) – Default is 1.
Returns:

UFL Form.

fenics_shells.reissner_mindlin.forms.psi_T(gamma, **kwargs)

Returns transverse shear energy density calculated using:

Isotropic case: .. math:

\psi_T(\gamma, \ldots) = \frac{E \kappa t}{4(1 + \nu))}\gamma**2
Parameters:
  • gamma – Shear strain, typically UFL form with shape (2,).
  • **kwargs – Isotropic case: E: Young’s modulus, Constant or Expression. nu: Poisson’s ratio, Constant or Expression. t: Thickness, Constant or Expression. kappa: Shear correction factor, Constant or Expression.
Returns:

UFL expression of transverse shear stress vector.

fenics_shells.reissner_mindlin.function_spaces module

fenics_shells.reissner_mindlin.function_spaces.DuranLibermanSpace(mesh)

A helper function which returns a FiniteElement for the simulation of the out-of-plane Reissner-Mindlin problem without shear-locking based on the ideas of Duran and Liberman’s paper:

R. Durán and E. Liberman, “On mixed finite element methods for the Reissner-Mindlin plate model” Math. Comp., vol. 58, no. 198, pp. 561–573, 1992.

Parameters:mesh (dolfin.Mesh) – a mesh of geometric dimension 2.
Returns:fenics_shells.ProjectedFunctionSpace.
fenics_shells.reissner_mindlin.function_spaces.MITC7Space(mesh, space_type='multipliers')

Warning: Currently not working due to regressions in FFC/FIAT.

A helper function which returns a FunctionSpace for the simulation of the out-of-plane Reissner-Mindlin problem without shear-locking based on the ideas of Brezzi, Bathe and Fortin’s paper:

“Mixed-interpolated elements for Reissner–Mindlin plates” Int. J. Numer. Meth. Engng., vol. 28, no. 8, pp. 1787–1801, Aug. 1989. http://dx.doi.org/10.1002/nme.1620280806

In the case that space_type is “multipliers”, a dolfin.FunctionSpace will be returned with an additional Lagrange multiplier field to tie the shear strains computer with the primal variables (transverse displacement and rotations) to the independent shear strain variable.

In the case that space_type is “primal”, a dolfin.FunctionSpace will be returned with just the primal (transverse displacement and rotation) variables. Of course, any Reissner-Mindlin problem constructed with this approach will be prone to shear-locking.

Parameters:
  • mesh (dolfin.Mesh) – a mesh of geometric dimension 2.
  • space_type (Optional[str]) – Can be “primal”, or
  • Default is "multipliers". ("multipliers".) –
Returns:

dolfin.FunctionSpace.

Module contents

Overview

This package contains routines for simulating thin through to moderately thick plate structures using the Reissner-Mindlin theory.

Further reading