Region#

class felupe.Region(mesh, element, quadrature, grad=True)[source]#

A numeric region as a combination of a mesh, an element and a numeric integration scheme (quadrature). The gradients of the element shape functions are evaluated at all integration points of each cell in the region if the optional gradient argument is True.

\[ \begin{align}\begin{aligned}\frac{\partial X^I}{\partial r^J} &= X_a^I \frac{\partial h_a}{\partial r^J}\\\frac{\partial h_a}{\partial X^J} &= \frac{\partial h_a}{\partial r^I} \frac{\partial r^I}{\partial X^J}\\dV &= \det\left(\frac{\partial X^I}{\partial r^J}\right) w\end{aligned}\end{align} \]
Parameters:
  • mesh (Mesh) – A mesh with points and cells.

  • element (Element) – The finite element formulation to be applied on the cells.

  • quadrature (Quadrature) – An element-compatible numeric integration scheme with points and weights.

  • grad (bool, optional) – A flag to invoke gradient evaluation (default is True).

mesh#

A mesh with points and cells.

Type:

Mesh

element#

The finite element formulation to be applied on the cells.

Type:

Finite element

quadrature#

An element-compatible numeric integration scheme with points and weights.

Type:

Quadrature scheme

h#

Element shape function array h_aq of shape function a evaluated at quadrature point q.

Type:

ndarray

dhdr#

Partial derivative of element shape function array dhdr_aJq with shape function a w.r.t. natural element coordinate J evaluated at quadrature point q for every cell c (geometric gradient or Jacobian transformation between X and r).

Type:

ndarray

dXdr#

Geometric gradient dXdr_IJqc as partial derivative of undeformed coordinate I w.r.t. natural element coordinate J evaluated at quadrature point q for every cell c (geometric gradient or Jacobian transformation between X and r).

Type:

ndarray

drdX#

Inverse of dXdr.

Type:

ndarray

dV#

Numeric Differential volume element as product of determinant of geometric gradient dV_qc = det(dXdr)_qc w_q and quadrature weight w_q, evaluated at quadrature point q for every cell c.

Type:

ndarray

dhdX#

Partial derivative of element shape functions dhdX_aJqc of shape function a w.r.t. undeformed coordinate J evaluated at quadrature point q for every cell c.

Type:

ndarray

Examples

>>> from felupe import Cube, Hexahedron, GaussLegendre, Region
>>> mesh = Cube(n=3)
>>> element = Hexahedron()
>>> quadrature = GaussLegendre(order=1, dim=3)
>>> region = Region(mesh, element, quadrature)
>>> region
<felupe Region object>
  Element formulation: Hexahedron
  Quadrature rule: GaussLegendre
  Gradient evaluated: True

Cell-volumes may be obtained by a sum of the differential volumes located at the quadrature points.

>>> region.dV.sum(axis=0)
array([0.125, 0.125, 0.125, 0.125, 0.125, 0.125, 0.125, 0.125])

The partial derivative of the first element shape function w.r.t. the undeformed coordinates evaluated at the second integration point of the last element of the region:

>>> region.dhdX[0, :, 1, -1]
array([-1.24401694, -0.33333333, -0.33333333])
class felupe.RegionBoundary(mesh, element, quadrature, grad=True, only_surface=True, mask=None, ensure_3d=False)[source]#

A numeric boundary-region as a combination of a mesh, an element and a numeric integration scheme (quadrature). The gradients of the element shape functions are evaluated at all integration points of each cell in the region if the optional gradient argument is True.

\[ \begin{align}\begin{aligned}\frac{\partial X^I}{\partial r^J} &= X_a^I \frac{\partial h_a}{\partial r^J}\\\frac{\partial h_a}{\partial X^J} &= \frac{\partial h_a}{\partial r^I} \frac{\partial r^I}{\partial X^J}\\dV &= \det\left(\frac{\partial X^I}{\partial r^J}\right) w\end{aligned}\end{align} \]
Parameters:
  • mesh (Mesh) – A mesh with points and cells.

  • element (Element) – The finite element formulation to be applied on the cells.

  • quadrature (Quadrature) – An element-compatible numeric integration scheme with points and weights.

  • grad (bool, optional) – A flag to invoke gradient evaluation (default is True).

  • only_surface (bool, optional) – A flag to use only the enclosing outline of the region (default is True).

  • mask (ndarray or None, optional) – A boolean array to select a specific set of points (default is None).

  • ensure_3d (bool, optional) – A flag to enforce 3d area normal vectors.

mesh#

A mesh with points and cells.

Type:

Mesh

element#

The finite element formulation to be applied on the cells.

Type:

Finite element

quadrature#

An element-compatible numeric integration scheme with points and weights.

Type:

Quadrature scheme

h#

Element shape function array h_aq of shape function a evaluated at quadrature point q.

Type:

ndarray

dhdr#

Partial derivative of element shape function array dhdr_aJq with shape function a w.r.t. natural element coordinate J evaluated at quadrature point q for every cell c (geometric gradient or Jacobian transformation between X and r).

Type:

ndarray

dXdr#

Geometric gradient dXdr_IJqc as partial derivative of undeformed coordinate I w.r.t. natural element coordinate J evaluated at quadrature point q for every cell c (geometric gradient or Jacobian transformation between X and r).

Type:

ndarray

drdX#

Inverse of dXdr.

Type:

ndarray

dA#

Numeric Differential area vectors.

Type:

ndarray

normals#

Area unit normal vectors.

Type:

ndarray

dV#

Numeric Differential volume element as norm of Differential area vectors.

Type:

ndarray

dhdX#

Partial derivative of element shape functions dhdX_aJqc of shape function a w.r.t. undeformed coordinate J evaluated at quadrature point q for every cell c.

Type:

ndarray

mesh_faces()[source]#

Return a Mesh with face-cells on the selected boundary.

class felupe.RegionQuad(mesh, quadrature=<felupe.quadrature._gausslegendre.GaussLegendre object>, grad=True)[source]#

Bases: Region

A region with a quad element.

class felupe.RegionQuadBoundary(mesh, quadrature=<felupe.quadrature._gausslegendre.GaussLegendreBoundary object>, grad=True, only_surface=True, mask=None, ensure_3d=False)[source]#

Bases: RegionBoundary

A region with a quad element.

class felupe.RegionHexahedronBoundary(mesh, quadrature=<felupe.quadrature._gausslegendre.GaussLegendreBoundary object>, grad=True, only_surface=True, mask=None)[source]#

Bases: RegionBoundary

A region with a hexahedron element.

class felupe.RegionHexahedron(mesh, quadrature=<felupe.quadrature._gausslegendre.GaussLegendre object>, grad=True)[source]#

Bases: Region

A region with a hexahedron element.

class felupe.RegionTriangle(mesh, quadrature=<felupe.quadrature._triangle.Triangle object>, grad=True)[source]#

Bases: Region

A region with a triangle element.

class felupe.RegionTetra(mesh, quadrature=<felupe.quadrature._tetra.Tetrahedron object>, grad=True)[source]#

Bases: Region

A region with a tetra element.

class felupe.RegionConstantQuad(mesh, quadrature=<felupe.quadrature._gausslegendre.GaussLegendre object>, grad=False)[source]#

Bases: Region

A region with a constant quad element.

class felupe.RegionConstantHexahedron(mesh, quadrature=<felupe.quadrature._gausslegendre.GaussLegendre object>, grad=False)[source]#

Bases: Region

A region with a constant hexahedron element.

class felupe.RegionQuadraticHexahedron(mesh, quadrature=<felupe.quadrature._gausslegendre.GaussLegendre object>, grad=True)[source]#

Bases: Region

A region with a (serendipity) quadratic hexahedron element.

class felupe.RegionTriQuadraticHexahedron(mesh, quadrature=<felupe.quadrature._gausslegendre.GaussLegendre object>, grad=True)[source]#

Bases: Region

A region with a tri-quadratic (lagrange) hexahedron element.

class felupe.RegionQuadraticTriangle(mesh, quadrature=<felupe.quadrature._triangle.Triangle object>, grad=True)[source]#

Bases: Region

A region with a quadratic triangle element.

class felupe.RegionQuadraticTetra(mesh, quadrature=<felupe.quadrature._tetra.Tetrahedron object>, grad=True)[source]#

Bases: Region

A region with a quadratic tetra element.

class felupe.RegionTriangleMINI(mesh, quadrature=<felupe.quadrature._triangle.Triangle object>, grad=True)[source]#

Bases: Region

A region with a triangle-MINI element.

class felupe.RegionTetraMINI(mesh, quadrature=<felupe.quadrature._tetra.Tetrahedron object>, grad=True)[source]#

Bases: Region

A region with a tetra-MINI element.

class felupe.RegionLagrange(mesh, order, dim, quadrature=None, grad=True)[source]#

Bases: Region

A region with an arbitrary order lagrange element.