Module shapes

Source
Expand description

Interpolation functions and derivatives for geometric shapes (elements)

§Definitions

Here, we consider the following definitions:

  • space_ndim – is the number of dimensions of the space under study (2 or 3)
  • geo_ndim – is the number of dimensions of the geometry element (shape). For instance, a line in the 2D space has geo_ndim = 1 and space_ndim = 2. Another example: a triangle in the 3D space has geo_ndim = 2 and space_ndim = 3.
  • local – refers to a numbering scheme for the nodes of the shape (or element)
  • global – refers to a numbering scheme applied for the whole mesh
  • spatial (real) space – is the “real” space mapped by the x₀,x₁,x₂ coordinates (see figure below)
  • reference (natural) space – is the “virtual” space mapped by the ξ₀,ξ₁,ξ₂ coordinates (see figure below)

Note: In the code, we use (r,s,t) as (ξ₀,ξ₁,ξ₂) and ksi to represent the vector ξ.

Real and reference space mapping

We also consider the following counting variables:

  • nnode – (local) number of points (aka nodes) that define the shape/element.
  • npoint – (global) number of points in the whole mesh; not used in this module but important to remember
  • nedge – number of edges on the shape (2D or 3D)
  • nface – number of faces on the shape (3D only)
  • edge_nnode – number of points/nodes that define the edge
  • face_nnode – number of points/nodes that define the face
  • face_nedge – number of edges on the face

Geometry cases regarding the number of dimensions (geo vs space)

  1. Case CABLEgeo_ndim = 1 and space_ndim = 2 or 3; e.g., line in 2D or 3D (cables and rods)
  2. Case SHELLgeo_ndim = 2 and space_ndim = 3; e.g. Tri or Qua in 3D (shells and surfaces)
  3. Case SOLIDgeo_ndim = space_ndim; e.g., Tri and Qua in 2D or Tet and Hex in 3D
geo_ndimspace_ndim = 2space_ndim = 3
1CABLECABLE
2SOLIDSHELL
3impossibleSOLID

§Isoparametric formulation

The isoparametric formulation establishes that we can calculate the coordinates x(ξ) within the shape/element from the shape functions Nᵐ(ξ) and the coordinates at each node by using the formula:

→ →         →  →
x(ξ) = Σ Nᵐ(ξ) xᵐ
       m

where x is the (space_ndim) vector of real coordinates, ξ is the (geo_ndim) vector of reference coordinates, Nᵐ are the (nnode) interpolation functions, and xᵐ are the (nnode) coordinates of each m-node of the geometric shape.

Given an (nnode,space_ndim) matrix of coordinates X, we can calculate the (space_ndim) vector of coordinates x by means of

x = Xᵀ ⋅ N

where N is an (nnode) vector formed with all Nᵐ values.

§Derivatives on the reference space and gradients on the real space

Here, we consider two cases:

  • General case with geo_ndim = space_ndim; and
  • Line in multi-dimensions with geo_ndim = 1 and space_ndim > 1.

§SOLID case with geo_ndim = space_ndim

If SOLID (geo_ndim = space_ndim = 2 or 3), we define the Jacobian tensor as

        →
  →    dx     →    →
J(ξ) = —— = Σ xᵐ ⊗ Lᵐ
        →   m
       dξ

where

            →
→  →    dNᵐ(ξ)
Lᵐ(ξ) = ——————
           →
          dξ

are the derivatives of each interpolation function Nᵐ with respect to the reference coordinate. Lᵐ are (geo_ndim) vectors and can be organized in an (nnode,geo_ndim) matrix L of local derivatives.

We can write the Jacobian in matrix (space_ndim,geo_ndim) notation as follows

J = Xᵀ · L

where X is the (nnode,space_ndim) matrix of coordinates and L is the (nnode,geo_ndim) matrix.

Next, we define the gradient of interpolation functions (i.e., derivatives of interpolation functions with respect to real coordinates) by means of

            →
→  →    dNᵐ(ξ)
Bᵐ(ξ) = ——————
           →
          dx

which can be organized in an (nnode,space_ndim) matrix B.

The inverse Jacobian allows us to determine the gradient vectors B as follows

→       →  →        →
Bᵐ(ξ) = Lᵐ(ξ) · J⁻¹(ξ)

Or, in matrix notation,

B = L · J⁻¹

where B is an (nnode,space_ndim) matrix.

§SHELL case with geo_ndim = 2 and space_ndim = 3

In this case, the Jacobian matrix is (3,2) and can also be computed by the following matrix multiplication

       dx
J(ξ) = ——
       dξ

Or, in matrix notation,

J = Jshell = Xᵀ · L

However, the inverse Jacobian and gradients are not available in this case.

§CABLE case with geo_ndim = 1 and space_ndim = 2 or 3

In this case, the Jacobian equals the (space_ndim,1) base vector g₁ which is tangent to the line element, i.e.,

                          →
→    →      →    →  →    dx
J := Jcable(ξ) = g₁(ξ) = ——
                         dξ
matrix notation: Jcable = Xᵀ · L

We also consider a parametric coordinate which varies from 0 to ℓ_max (the length of the line) according to

               ℓ_max
ℓ(ξ) = (1 + ξ) —————
                 2

       2 · ℓ
ξ(ℓ) = ————— - 1
       ℓ_max

Note that:

0 ≤ ℓ ≤ ℓ_max

-1 ≤ ξ ≤ +1

§Normal vectors

§CABLE case with geo_ndim = 1 and space_ndim = 2 or 3

Base vector tangent to the line:

         →
        dx
g₁(ξ) = —— = Xᵀ · L = first_column(J)
        dξ

Normal vector:

→   →    →
n = e₃ × g₁

  →       →
||n|| = ||g₁||

Thus

       →           →
dℓ = ||g₁|| dξ = ||n|| dξ

For a straight line (segment):

  →
||n|| = Δℓ / Δξ = L / 2

because all GeoClass::Lin have Δξ = 2.

§SHELL case with geo_ndim = 2 and space_ndim = 3

Base vectors tangent to the surface:

         →
→  →    dx
g₁(ξ) = ——— = first_column(J)
        dξ₁

         →
→  →    dx
g₂(ξ) = ——— = second_column(J)
        dξ₂

Normal vector:

→   →    →
n = g₁ × g₂

Thus

        →
dA := ||n|| dξ₁ dξ₂

For flat quadrilateral faces with sides perpendicular one with another

  →
||n|| = A / (Δξ₁ Δξ₂) = A / 4

because all GeoClass::Qua have Δξᵢ = 2.

§GeoKind and GeoClass

GeoKind is perhaps the most important enum in this module. It defines the actual shape used in finite element analyses.

Below are some example of shapes, classified according to GeoClass. The numbers are the local numbers (nodes).

§Lines – Lin

lin_cells

§Triangles – Tri

tri_cells

§Quadrilaterals – Qua

qua_cells

§Tetrahedra – Tet

tet_cells

§Hexahedra – Hex

hex_cells

Structs§

Hex8
Defines a hexahedron with 8 nodes (bilinear faces)
Hex20
Defines a hexahedron with 20 nodes (quadratic faces)
Hex32
Defines a hexahedron with 32 nodes (cubic faces)
Lin2
Defines a line (segment) with 2 nodes (linear functions)
Lin3
Defines a line (segment) with 3 nodes (quadratic functions)
Lin4
Defines a line (segment) with 4 nodes (cubic functions)
Lin5
Defines a line (segment) with 5 nodes (quartic functions)
Qua4
Defines a quadrilateral with 4 nodes (linear edges)
Qua8
Defines a quadrilateral with 8 nodes (quadratic edges)
Qua9
Defines a quadrilateral with 9 nodes (quadratic edges; interior node)
Qua12
Defines a quadrilateral with 12 nodes (cubic edges)
Qua16
Defines a quadrilateral with 16 nodes (cubic edges; interior nodes)
Qua17
Defines a quadrilateral with 17 nodes (quartic edge; interior node)
Scratchpad
Holds the variables used by the functions in the sub-module named op
Tet4
Defines a tetrahedron with 4 nodes (linear faces)
Tet10
Defines a tetrahedron with 10 nodes (quadratic faces)
Tet20
Defines a tetrahedron with 20 nodes (cubic faces)
Tri3
Defines a triangle with 3 nodes (linear edges)
Tri6
Defines a triangle with 6 nodes (quadratic edges)
Tri10
Defines a triangle with 10 nodes (cubic edges; interior node)
Tri15
Defines a triangle with 15 nodes (quartic edges; interior nodes)

Enums§

GeoCase
Defines the geometry case regarding the number of dimensions (geo vs space) (Cable, Shell, Solid)
GeoClass
Defines the class of the geometric shape (Lin, Tri, Qua, Tet, Hex)
GeoKind
Defines the kind of the geometric shape (Lin2, … Tri3, …, Qua4, … Tet4, …, Hex8, …)

Constants§

DET_JAC_NOT_AVAILABLE
Indicates that the determinant of the Jacobian is not available (e.g., Shells)
HEX_EDGE_TO_FACE
Converts edge-local-index to face-local-indices of a Hexahedron

Functions§

geo_case
Returns the geometry case given the geo and space dimensions

Type Aliases§

FnDeriv
Defines a function to calculate the derivatives of interpolation functions
FnInterp
Defines a function to calculate interpolation functions