fenics_shells.fem package

Submodules

fenics_shells.fem.CDG module

fenics_shells.fem.CDG.cdg_energy(theta, M, stabilization, mesh, bcs_theta=None, dS=<Mock id='140702363087184'>)

Return the continuous/discontinuous terms for a fourth-order plate model.

\[\pi_{cdg} = - \partial_n w \cdot M_{n}(w) + \frac{1}{2} \frac{\alpha}{|e|} |\partial_n w |^2\]
Parameters:
  • theta – Rotations, UFL or DOLFIN Function of rank (2,) (vector).
  • M – UFL form of bending moment tensor of rank (2,2) (tensor).
  • stabilization – a constant or ulf expression providing the stabilization parameter of the continuous/discontinuous formulation. This should be an eximation of the norm of the bending stiffness
  • mesh – DOLFIN mesh.
  • bcs_theta (Optional) – list of dolfin.DirichletBC for the rotations theta. Defaults to None.
  • dS – (Optional). Measure on interior facets. Defaults to dolfin.dS.
Returns:

a dolfin.Form associated with the continuous/discontinuous formulation.

The Kirchhoff-Love plate model is a fourth-order PDE, giving rise to a weak form with solution in Sobolev space \(H^2(\Omega)\). Because FEniCS does not currently include support for \(H^2(\Omega)\) conforming elements we implement a hybrid continuous/discontinuous approach, allowing the use of Lagrangian elements with reduced regularity requirements.

Description can be found in the paper:
G. Engel, K. Garikipati, T. J. R. Hughes, M. G. Larson, L. Mazzei and R. L. Taylor, “Continuous/discontinuous finite element approximations of fourth-order elliptic problems in structural and continuum mechanics with applications to thin beams and plates, and strain gradient elasticity” Comput. Method. Appl. M., vol. 191, no. 34, pp. 3669-3750, 2002.
fenics_shells.fem.CDG.cdg_stabilization(E, t)

Returns the stabilization parameter as the norm of the bending stiffness matrix.

Parameters:
  • E – Young’s modulus, Constant or Expression.
  • t – Thickness, Constant or Expression.
Returns:

a dolfin.Coefficient providing the stabilization parameter of the continuous/discontinuous formulation.

fenics_shells.fem.assembling module

fenics_shells.fem.assembling.assemble(*args, **kwargs)

Pass-through for dolfin.assemble and fenics_shells.projected_assemble.

If the first argument is an instance of ProjectedFunctionSpace it will call fenics_shells.projected_assemble, otherwise it will pass through to dolfin.assemble.

fenics_shells.fem.assembling.projected_assemble(U_P, a, L, bcs=None, A=None, b=None, is_interpolation=False, a_is_symmetric=False, form_compiler_parameters=None, add_values=False, finalize_tensor=True, keep_diagonal=False, backend=None)

fenics_shells.fem.solving module

fenics_shells.fem.solving.reconstruct_full_space(u_f, u_p, a, L, is_interpolation=False, a_is_symmetric=False, form_compiler_parameters=None)

Given a Function on a projected space \(u_p \in U_P\) and a function in the full space \(u_f \in U_F: such that :math:`U_P \subset U_F\), reconstruct the variable u_f on the full space via direct copy of the matching subfunctions shared between U_F and U_P and the local solution of the original problem a == L.

Parameters:
  • u_f – DOLFIN Function on FullFunctionSpace.
  • u_p – DOLFIN Function on ProjectedFunctionSpace.
  • a – Bilinear form.
  • L – Bilinear form.
Returns:

DOLFIN Function on FullFunctionSpace.

Return type:

u_f

Module contents

fenics_shells.fem.projected_assemble(U_P, a, L, bcs=None, A=None, b=None, is_interpolation=False, a_is_symmetric=False, form_compiler_parameters=None, add_values=False, finalize_tensor=True, keep_diagonal=False, backend=None)
fenics_shells.fem.assemble(*args, **kwargs)

Pass-through for dolfin.assemble and fenics_shells.projected_assemble.

If the first argument is an instance of ProjectedFunctionSpace it will call fenics_shells.projected_assemble, otherwise it will pass through to dolfin.assemble.

fenics_shells.fem.reconstruct_full_space(u_f, u_p, a, L, is_interpolation=False, a_is_symmetric=False, form_compiler_parameters=None)

Given a Function on a projected space \(u_p \in U_P\) and a function in the full space \(u_f \in U_F: such that :math:`U_P \subset U_F\), reconstruct the variable u_f on the full space via direct copy of the matching subfunctions shared between U_F and U_P and the local solution of the original problem a == L.

Parameters:
  • u_f – DOLFIN Function on FullFunctionSpace.
  • u_p – DOLFIN Function on ProjectedFunctionSpace.
  • a – Bilinear form.
  • L – Bilinear form.
Returns:

DOLFIN Function on FullFunctionSpace.

Return type:

u_f