Element types

Generic classes

class stagpyviz.Element

Base class for isoparametric elements. This class defines the interface for evaluating shape functions, their derivatives, and Jacobian-related quantities.

This class is the parent class for some specific element types.

Note

In the following documentation, the term reference refers to the reference element (e.g., the unit triangle or unit square), while physical refers to the actual element in the physical domain defined by its node coordinates.

evaluate_Jacobian(GNi: ndarray, xe: ndarray) ndarray

Compute the Jacobian matrix \(\boldsymbol J\) for a single element or multiple elements such that:

\[J_{ij} = \sum_k \frac{\partial N_k}{\partial \xi_j} x_{k,i}\]

where \(\partial_{\xi_j} N_k\) are the derivatives of the shape functions with respect to the reference coordinates, and \(x_{k,i}\) are the physical coordinates of the node \(k\).

Parameters:
  • GNi (numpy.ndarray) – Array of shape (nnodes, reference dim) containing the derivatives of the shape functions with respect to the reference coordinates.

  • xe (numpy.ndarray) – Array of shape (nnodes, physical dim) for a single element, or (n_cells, nnodes, physical dim) for multiple elements, containing the physical coordinates of the nodes.

Returns:

The Jacobian matrix \(\boldsymbol J\) for the element(s). For a single element, array of the shape (reference dim, physical dim). For multiple elements, array of the shape (n_cells, reference dim, physical dim).

Return type:

numpy.ndarray

evaluate_dNidx(invJ: ndarray, GNi: ndarray) ndarray

Compute the physical derivatives of the shape functions with respect to the physical coordinates such that:

\[\frac{\partial N_k}{\partial x_i} = \sum_j J^{-1}_{ji} \frac{\partial N_k}{\partial \xi_j}\]

where \(J^{-1}_{ji}\) are the components of the inverse Jacobian matrix (transposed), and \(\partial_{\xi_j} N_k\) are the derivatives of the shape functions with respect to the reference coordinates.

Parameters:
  • invJ (numpy.ndarray) – Array of shape (reference dim, physical dim) for a single element, or (n_cells, reference dim, physical dim) for multiple elements, containing the inverse Jacobian matrix (transposed).

  • GNi (numpy.ndarray) – Array of shape (nnodes, dim) containing the derivatives of the shape functions with respect to the reference coordinates.

Returns:

The physical derivatives of the shape functions with respect to the physical coordinates. For a single element, array of the shape (nnodes, physical dim). For multiple elements, array of the shape (n_cells, nnodes, physical dim).

Return type:

numpy.ndarray

evaluate_element_centroid(xe: ndarray) ndarray

Compute the element centroid in physical coordinates such that:

\[c_{i} = \sum_k N_k x_{k,i}\]

where \(N_k\) are the shape functions evaluated at the element centroid, and \(x_{k,i}\) are the physical coordinates of the node \(k\).

Parameters:

xe (numpy.ndarray) – Array of shape (nnodes, physical dim) for a single element, or (n_cells, nnodes, physical dim) for multiple elements, containing the physical coordinates of the nodes.

Returns:

The element centroid in physical coordinates. For a single element, an array of shape (physical dim,). For multiple elements, an array of shape (n_cells, physical dim).

Return type:

numpy.ndarray

evaluate_volume(xe: ndarray, rule: int = 1) ndarray

Evaluate the volume (or area in 2D) of the element(s) using numerical quadrature such that:

\[V_e = \int_{\Omega_e} 1 \, dV \approx \sum_q w_q |\det(\boldsymbol J_q)|\]

where \(w_q\) are the quadrature weights, \(\det(\boldsymbol J_q)\) are the determinants of the Jacobian matrices at the quadrature points, and the sum is over all quadrature points.

Parameters:
  • xe (numpy.ndarray) – Array of shape (n_cells, nnodes, physical dim) containing the physical coordinates of the nodes.

  • rule (int) – The quadrature rule to use. See elements specific documentation.

Returns:

The volume (or area in 2D) of the element(s). Array of shape (n_cells,).

Return type:

numpy.ndarray

integrate_field(xe: ndarray, field_e: ndarray, rule: int = 1) ndarray

Evaluate the integral of a field defined at the nodes of the element(s) using numerical quadrature such that:

\[I_e = \int_{\Omega_e} f(\mathbf x) \, dV \approx \sum_q w_q |\det(\boldsymbol J_q)| \left( \sum_k N_k(\boldsymbol \xi_q) f_k \right)\]

where \(w_q\) are the quadrature weights, \(\det(\boldsymbol J_q)\) are the determinants of the Jacobian matrices at the quadrature points, \(N_k(\boldsymbol \xi_q)\) are the shape functions evaluated at the quadrature points and \(f_k\) are the values of the field at the nodes.

Parameters:
  • xe (numpy.ndarray) – Array of shape (n_cells, nnodes, physical dim) containing the physical coordinates of the nodes.

  • field_e (numpy.ndarray) – Array of shape (n_cells, nnodes, ...) containing the values of the field at the nodes. The last dimensions can be of any shape, and will be preserved in the output.

  • rule (int) – The quadrature rule to use. See elements specific documentation.

Returns:

The integral of the field over the element(s). Array of shape (n_cells, ...).

Return type:

numpy.ndarray

class stagpyviz.Element2D

Class for 2D elements (e.g., triangles, quadrilaterals) representing domains in 2D space, \(\mathbb R^2\).

evaluate_detJ(J: ndarray) ndarray

Compute the determinant of the Jacobian matrix \(\boldsymbol J\) for a single element or multiple elements such that:

\[\det(\boldsymbol J) = J_{00} J_{11} - J_{01} J_{10}\]
Parameters:

J (numpy.ndarray) – Array of shape (2, 2) for a single element, or (n_cells, 2, 2) for multiple elements, containing the Jacobian matrix.

Returns:

The determinant of the Jacobian matrix \(\det(\boldsymbol J)\) for the element(s). For a single element, a float. For multiple elements, an array of shape (n_cells,).

Return type:

float or numpy.ndarray

evaluate_invJ(J: ndarray, detJ: float | ndarray) ndarray

Compute the inverse of the Jacobian matrix \(\boldsymbol J\) for a single element or multiple elements such that:

\[\begin{split}\boldsymbol J^{-1} = \frac{1}{\det(\boldsymbol J)} \begin{bmatrix} J_{11} & -J_{01} \\ -J_{10} & J_{00} \end{bmatrix}\end{split}\]
Parameters:
  • J (numpy.ndarray) – Array of shape (2, 2) for a single element, or (n_cells, 2, 2) for multiple elements, containing the Jacobian matrix.

  • detJ (float or numpy.ndarray) – The determinant of the Jacobian matrix for the element(s). For a single element, a float. For multiple elements, an array of shape (n_cells,).

Returns:

The inverse of the Jacobian matrix \(\boldsymbol J^{-1}\) for the element(s). For a single element, an array of shape (2, 2). For multiple elements, an array of shape (n_cells, 2, 2).

Return type:

numpy.ndarray

class stagpyviz.Element3D

Class for 3D elements (e.g., tetrahedra, hexahedra) representing domains in 3D space, \(\mathbb R^3\).

evaluate_detJ(J: ndarray) float | ndarray

Compute the determinant of the Jacobian matrix \(\boldsymbol J\) for a single element or multiple elements such that:

\[\det(\boldsymbol J) = J_{00}(J_{11}J_{22} - J_{12}J_{21}) - J_{01}(J_{10}J_{22} - J_{12}J_{20}) + J_{02}(J_{10}J_{21} - J_{11}J_{20})\]
Parameters:

J (numpy.ndarray) – Array of shape (..., 3, 3) containing the Jacobian matrix, where the last two dimensions are the matrix indices.

Returns:

The determinant of the Jacobian matrix \(\det(\boldsymbol J)\) for the element(s). For a single element, a float. For multiple elements, an array with shape J.shape[:-2].

Return type:

float or numpy.ndarray

evaluate_invJ(J: ndarray, detJ: float | ndarray) ndarray

Compute the inverse of the Jacobian matrix \(\boldsymbol J\) for a single element or multiple elements.

Parameters:
  • J (numpy.ndarray) – Array of shape (3, 3) for a single element, or (n_cells, 3, 3) for multiple elements, containing the Jacobian matrix.

  • detJ (float or numpy.ndarray) – The determinant of the Jacobian matrix for the element(s). For a single element, a float. For multiple elements, an array of shape (n_cells,).

Returns:

The inverse of the Jacobian matrix \(\boldsymbol J^{-1}\) for the element(s). For a single element, an array of shape (3, 3). For multiple elements, an array of shape (n_cells, 3, 3).

Return type:

numpy.ndarray

class stagpyviz.SurfaceElement

Specialized class for 2D surface elements (triangles, quadrilaterals) embedded in 3D space, \(\mathbb R^3\).

evaluate_dNidx(J: ndarray, GNi: ndarray) ndarray

Compute the physical derivatives of the shape functions with respect to the physical coordinates for a surface element such that:

\[\frac{\partial N_k}{\partial x_i} = \sum_l \left( \sum_m J_{im} G^{-1}_{ml} \right) \frac{\partial N_k}{\partial \xi_l}\]

where \(J_{im}\) are the components of the Jacobian matrix, \(G^{-1}_{ml}\) are the components of the inverse of the metric tensor, and \(\partial_{\xi_l} N_k\) are the derivatives of the shape functions with respect to the reference coordinates.

Parameters:
  • J (numpy.ndarray) – Array of shape (..., 3, 2) containing the Jacobian matrix. The ellipsis can represent multiple elements and/or multiple evaluation points, but the last two dimensions must be (3, 2).

  • GNi (numpy.ndarray) – Array of shape (..., nnodes, 2) containing the derivatives of the shape functions with respect to the reference coordinates. The ellipsis can represent multiple elements and/or multiple evaluation points, but the last two dimensions must be (nnodes, 2).

Returns:

The physical derivatives of the shape functions with respect to the physical coordinates. Shape (..., nnodes, 3). The ellipsis can represent multiple elements and/or multiple evaluation points depending on the input shapes.

Return type:

numpy.ndarray

evaluate_detG(J: ndarray) float | ndarray

Compute the determinant of the metric tensor \(\boldsymbol G\) for a single element or multiple elements such that:

\[\det(\boldsymbol G) = \mathbf n \cdot \mathbf n\]

where \(\mathbf n\) is the non-unit normal vector to the surface element computed using the method normal_vector_nonu.

Parameters:

J (numpy.ndarray) – Array of shape (..., 3, 2) containing the Jacobian matrix. The ellipsis can represent multiple elements and/or multiple evaluation points, but the last two dimensions must be (3, 2).

Returns:

The determinant of the metric tensor \(\det(\boldsymbol G)\) for the element(s). Shape of the input ellipsis.

Return type:

float or numpy.ndarray

evaluate_detJ(J: ndarray) float | ndarray

Compute the determinant of the Jacobian matrix \(\boldsymbol J\) for a surface element such that:

\[\det(\boldsymbol J) = \sqrt{\det(\boldsymbol G)}\]

where \(\det(\boldsymbol G)\) is the determinant of the metric tensor computed using the method evaluate_detG.

Parameters:

J (numpy.ndarray) – Array of shape (..., 3, 2) containing the Jacobian matrix. The ellipsis can represent multiple elements and/or multiple evaluation points, but the last two dimensions must be (3, 2).

Returns:

The determinant of the Jacobian matrix for the element(s). Shape of the input ellipsis.

Return type:

float or numpy.ndarray

evaluate_metric_tensor(J: ndarray) ndarray

Compute the metric tensor \(\boldsymbol G\) for a surface element such that:

\[G_{ij} = \mathbf t_{\xi_i} \cdot \mathbf t_{\xi_j}\]

where \(\mathbf t_{\xi_i}\) and \(\mathbf t_{\xi_j}\) are the tangent vectors defined by the columns of the Jacobian matrix.

Parameters:

J (numpy.ndarray) – Array of shape (..., 3, 2) containing the Jacobian matrix. The ellipsis can represent multiple elements and/or multiple evaluation points, but the last two dimensions must be (3, 2).

Returns:

The metric tensor \(\boldsymbol G\) for the element(s). Array of shape (..., 2, 2). The ellipsis can represent multiple elements and/or multiple evaluation points depending on the input shape.

Return type:

numpy.ndarray

normal_vector(J: ndarray) ndarray

Compute the unit normal vector to the surface element. The non-unit normal vector is computed using the method normal_vector_nonu, and then normalized by its magnitude (which is the determinant of the Jacobian matrix) to obtain the unit normal vector.

Parameters:

J (numpy.ndarray) – Array of shape (..., 3, 2) containing the Jacobian matrix. The ellipsis can represent multiple elements and/or multiple evaluation points, but the last two dimensions must be (3, 2).

Returns:

The unit normal vector to the surface element. Array of shape (..., 3).

Return type:

numpy.ndarray

normal_vector_nonu(J: ndarray) ndarray

Compute the non-unit normal vector to the surface element defined by the Jacobian matrix \(\boldsymbol J\) for a single element or multiple elements. The normal vector \(\mathbf n\) is computed as the cross product of the two tangent vectors defined by the columns of the Jacobian matrix:

\[\begin{split}\begin{split} \boldsymbol J &= \begin{bmatrix} \mathbf t_{\xi} & \mathbf t_{\eta} \end{bmatrix} \in \mathbb R^{3 \times 2} \\ \mathbf n &= \mathbf t_{\xi} \times \mathbf t_{\eta} \end{split}\end{split}\]

where \(\mathbf t_{\xi}\) and \(\mathbf t_{\eta}\) are the tangent vectors defined by the columns of the Jacobian matrix, and \(\times\) denotes the cross product.

Parameters:

J (numpy.ndarray) – Array of shape (..., 3, 2) containing the Jacobian matrix. The ellipsis can represent multiple elements and/or multiple evaluation points, but the last two dimensions must be (3, 2).

Returns:

The non-unit normal vector to the surface element. Array of shape (..., 3).

Return type:

numpy.ndarray

Specific element types

class stagpyviz.P1_2D

Triangular linear \(\mathcal{P}_1\) parametric element in 2D, \(\mathbb R^2\). Nodes are ordered couterclockwise and the reference coordinates \(\xi, \eta\) range from 0 to 1 with the constraint \(\xi + \eta \leq 1\). Inherits from Element2D.

P1 triangle reference element
GNi_centroid() ndarray

Computes the shape function derivatives at the element centroid. The derivatives of the shape functions with respect to the reference coordinates are constant for a linear triangle, so this method simply returns the same values as evaluate_GNi.

Returns:

Array of shape (3, 2) containing the derivatives of the shape functions with respect to the reference coordinates at the element centroid.

Return type:

numpy.ndarray

Ni_centroid()

Compute the shape function values at the element centroid. For the chosen linear triangle, the shape functions at the centroid are all equal to 1/3.

Returns:

Array of shape (3,) containing the shape function values at the element centroid.

Return type:

numpy.ndarray

evaluate_GNi(xi: ndarray) ndarray

Compute the derivatives of the shape functions with respect to the reference coordinates. For a linear triangle, the derivatives are constant and given by:

\[\begin{split}\frac{\partial N_k}{\partial \xi_i} = \begin{bmatrix} -1 & -1 \\ 1 & 0 \\ 0 & 1 \end{bmatrix}\end{split}\]
Returns:

Array of shape (3, 2) containing the derivatives of the shape functions with respect to the reference coordinates.

Return type:

numpy.ndarray

evaluate_Ni(xi: ndarray)

Evaluate the shape functions at given reference coordinates \(\xi = (\xi, \eta)\). Shape functions for a linear triangle are:

\[\begin{split}\begin{split} N_0(\xi, \eta) &= 1 - \xi - \eta \\ N_1(\xi, \eta) &= \xi \\ N_2(\xi, \eta) &= \eta \end{split}\end{split}\]
Parameters:

xi (numpy.ndarray) – Array of shape (2,) containing the reference coordinates.

Returns:

Array of shape (3,) containing the shape function values at the given reference coordinates.

Return type:

numpy.ndarray

quadrature(nqp: int)

Return the quadrature weights and points for the specified number of quadrature points. Supported values for nqp are:

Parameters:

nqp (int) – Number of quadrature points. Supported values are 1 and 3.

Returns:

Tuple containing the quadrature weights and points.

Return type:

tuple(numpy.ndarray, numpy.ndarray)

quadrature_1pt()

1 point quadrature rule for a triangular element in 2D, which corresponds to evaluating the integrand at the centroid of the reference element.

\[\begin{split}w_0 = \frac{1}{2}, \quad \boldsymbol \xi_0 = \begin{bmatrix} 1/3 \\ 1/3 \end{bmatrix}\end{split}\]

with \(w_0\) the quadrature weight and \(\boldsymbol \xi_0\) the quadrature point.

Returns:

Tuple containing the quadrature weights and points. The weight is returned as a 1D array of shape (1,), and the points are returned as a 2D array of shape (1, 2).

Return type:

tuple(numpy.ndarray, numpy.ndarray)

quadrature_3pt()

3 point quadrature rule for a triangular element in 2D, which corresponds to evaluating the integrand at the three points defined by:

\[\begin{split}\begin{split} w_q &= \frac{1}{6} \\ \boldsymbol \xi_0 &= \begin{bmatrix} 1/6 \\ 1/6 \end{bmatrix} \quad \boldsymbol \xi_1 = \begin{bmatrix} 2/3 \\ 1/6 \end{bmatrix} \quad \boldsymbol \xi_2 = \begin{bmatrix} 1/6 \\ 2/3 \end{bmatrix} \end{split}\end{split}\]

with \(w_q\) the quadrature weights (all equal to 1/6 for this rule) and \(\boldsymbol \xi_q\) the quadrature points.

Returns:

Tuple containing the quadrature weights and points. The weights are returned as a 1D array of shape (3,), and the points are returned as a 2D array of shape (3, 2).

Return type:

tuple(numpy.ndarray, numpy.ndarray)

class stagpyviz.P1_2D_R3

Triangular linear \(\mathcal{P}_1\) parametric element in 3D, representing a surface element in \(\mathbb R^3\). Nodes are ordered couterclockwise and the reference coordinates \(\xi, \eta\) range from 0 to 1 with the constraint \(\xi + \eta \leq 1\). Inherits from SurfaceElement.

GNi_centroid() ndarray

Compute the shape function derivatives at the element centroid. See GNi_centroid for the shape function derivatives at the centroid.

Ni_centroid()

Compute the shape function values at the element centroid. See Ni_centroid for the shape function values at the centroid.

evaluate_GNi(xi: ndarray) ndarray

Compute the derivatives of the shape functions with respect to the reference coordinates. See evaluate_GNi for the shape function derivatives with respect to the reference coordinates.

evaluate_Ni(xi: ndarray)

Evaluate the shape functions at given reference coordinates \(\xi = (\xi, \eta)\). See evaluate_Ni for the shape function definitions.

quadrature(nqp: int)

Return the quadrature weights and points for the specified number of quadrature points. Supported values for nqp are:

Parameters:

nqp (int) – Number of quadrature points. Supported values are 1 and 3.

Returns:

Tuple containing the quadrature weights and points.

Return type:

tuple(numpy.ndarray, numpy.ndarray)

class stagpyviz.Q1_2D

Quadrilateral bilinear \(\mathcal{Q}_1\) parametric element in 2D, \(\mathbb R^2\). Reference coordinates \(\boldsymbol \xi = (\xi, \eta) \in [-1, 1]^2\). Inherits from Element2D.

Q1 quadrilateral reference element
GNi_centroid() ndarray

Evaluate the derivatives of the shape functions with respect to the reference coordinates at the element centroid. See evaluate_GNi for the shape function derivatives with respect to the reference coordinates.

Returns:

Array of shape (4, 2) containing the derivatives of the shape functions with respect to the reference coordinates at the element centroid.

Return type:

numpy.ndarray

Ni_centroid() ndarray

Return the shape function values at the centroid of the element. For the chosen bilinear quadrilateral, the shape functions at the centroid are all equal to 1/4.

Returns:

Array of shape (4,) containing the shape function values at the element centroid.

Return type:

numpy.ndarray

evaluate_GNi(xi: ndarray) ndarray

Evaluate the derivatives of the shape functions with respect to the reference coordinates at given reference coordinates \(\boldsymbol \xi = (\xi, \eta)\). The derivatives of the shape functions for a bilinear quadrilateral are given by:

\[\begin{split}\frac{\partial N_k}{\partial \boldsymbol \xi} = \frac{1}{4} \begin{bmatrix} -(1 - \eta) & -(1 - \xi) \\ (1 - \eta) & -(1 + \xi) \\ -(1 + \eta) & -(1 - \xi) \\ (1 + \eta) & -(1 + \xi) \end{bmatrix}\end{split}\]

where the rows correspond to the shape functions \(N_k\) and the columns correspond to the reference coordinates \((\xi, \eta)\).

Parameters:

xi (numpy.ndarray) – Array of shape (2,) containing the reference coordinates.

Returns:

Array of shape (4, 2) containing the derivatives of the shape functions with respect to the reference coordinates at the given reference coordinates.

Return type:

numpy.ndarray

evaluate_Ni(xi: ndarray) ndarray

Evaluate the shape functions at given reference coordinates \(\boldsymbol \xi = (\xi, \eta)\). The shape functions for a bilinear quadrilateral are given by:

\[\begin{split}\begin{split} N_0(\xi, \eta) &= \frac{1}{4} (1 - \xi) (1 - \eta) \\ N_1(\xi, \eta) &= \frac{1}{4} (1 + \xi) (1 - \eta) \\ N_2(\xi, \eta) &= \frac{1}{4} (1 - \xi) (1 + \eta) \\ N_3(\xi, \eta) &= \frac{1}{4} (1 + \xi) (1 + \eta) \end{split}\end{split}\]
Parameters:

xi (numpy.ndarray) – Array of shape (2,) containing the reference coordinates.

Returns:

Array of shape (4,) containing the shape function values at the given reference coordinates.

Return type:

numpy.ndarray

quadrature(nqp) tuple[ndarray, ndarray]

Return the quadrature weights and points for the specified number of quadrature points. Supported values for nqp are:

Parameters:

nqp (int) – Number of quadrature points. Supported values are 1 and 4.

Returns:

Tuple containing the quadrature weights and points.

Return type:

tuple(numpy.ndarray, numpy.ndarray)

quadrature_1pt() tuple[ndarray, ndarray]

1 point quadrature rule for a quadrilateral element in 2D, which corresponds to evaluating the integrand at the centroid of the reference element.

\[\begin{split}w_0 = 4, \quad \boldsymbol \xi_0 = \begin{bmatrix} 0 \\ 0 \end{bmatrix}\end{split}\]

with \(w_0\) the quadrature weight and \(\boldsymbol \xi_0\) the quadrature point.

Returns:

Tuple containing the quadrature weights and points. The weight is returned as a 1D array of shape (1,), and the points are returned as a 2D array of shape (1, 2).

Return type:

tuple(numpy.ndarray, numpy.ndarray)

quadrature_4pt() tuple[ndarray, ndarray]

\(2 \times 2\) point quadrature rule for a quadrilateral element in 2D.

\[\begin{split}w_q = 1, \quad \boldsymbol \xi_q = \begin{bmatrix} \pm 1/\sqrt{3}\\ \pm 1/\sqrt{3} \end{bmatrix}\end{split}\]

with \(w_q\) the quadrature weights (all equal to 1 for this rule) and \(\boldsymbol \xi_q\) the quadrature points.

Returns:

Tuple containing the quadrature weights and points. The weights are returned as a 1D array of shape (4,), and the points are returned as a 2D array of shape (4, 2).

Return type:

tuple(numpy.ndarray, numpy.ndarray)

class stagpyviz.Wedge3D

Wedge (prism) linear \(\mathcal{W}_1\) parametric element in 3D, \(\mathbb R^3\). Composed of two triangular faces and three quadrilateral faces. Inherits from Element3D.

Wedge reference element
GNi_centroid()

Evaluate the derivatives of the shape functions with respect to the reference coordinates at the element centroid. Calls the method evaluate_GNi with the reference coordinates of the centroid \((\xi, \eta, \zeta) = (1/3, 1/3, 0)\).

Returns:

Array of shape (6, 3) containing the derivatives of the shape functions with respect to the reference coordinates at the element centroid.

Return type:

numpy.ndarray

Ni_centroid()

Return the shape function values at the centroid of the element. For a linear wedge, the shape functions at the centroid are all equal to 1/6.

Returns:

Array of shape (6,) containing the shape function values at the element centroid.

Return type:

numpy.ndarray

evaluate_GNi(xi: ndarray)

Evaluate the derivatives of the shape functions with respect to the reference coordinates at given reference coordinates \(\boldsymbol \xi = (\xi, \eta, \zeta)\).

\[\begin{split}\frac{\partial N_k}{\partial \boldsymbol \xi} &= \frac{1}{2} \begin{bmatrix} -(1 - \zeta) & -(1 - \zeta) & -(1 - \xi - \eta) \\ 1 - \zeta & 0 & -\xi \\ 0 & 1 - \zeta & -\eta \\ -(1 + \zeta) & -(1 + \zeta) & 1 - \xi - \eta \\ 1 + \zeta & 0 & \xi \\ 0 & 1 + \zeta & \eta \\ \end{bmatrix} \\\end{split}\]

where the rows correspond to the shape functions \(N_k\) and the columns correspond to the reference coordinates \((\xi, \eta, \zeta)\).

Parameters:

xi (numpy.ndarray) – Array of shape (npoints, 3) containing the reference coordinates.

Returns:

Array of shape (npoints, 6, 3) containing the derivatives of the shape functions with respect to the reference coordinates at the given reference coordinates.

Return type:

numpy.ndarray

evaluate_Ni(xi: ndarray)

Evaluate the shape functions at given reference coordinates \(\boldsymbol \xi = (\xi, \eta, \zeta)\). The shape functions for a linear wedge are given by the product of the triangular shape functions in the \((\xi, \eta)\) plane and the linear shape functions in the \(\zeta\) direction:

\[\begin{split}\begin{split} N_0(\xi, \eta, \zeta) &= (1 - \xi - \eta) \cdot \frac{1 - \zeta}{2} \\ N_1(\xi, \eta, \zeta) &= \xi \cdot \frac{1 - \zeta}{2} \\ N_2(\xi, \eta, \zeta) &= \eta \cdot \frac{1 - \zeta}{2} \\ N_3(\xi, \eta, \zeta) &= (1 - \xi - \eta) \cdot \frac{1 + \zeta}{2} \\ N_4(\xi, \eta, \zeta) &= \xi \cdot \frac{1 + \zeta}{2} \\ N_5(\xi, \eta, \zeta) &= \eta \cdot \frac{1 + \zeta}{2} \end{split}\end{split}\]
Parameters:

xi (numpy.ndarray) – Array of shape (3,) containing the reference coordinates.

Returns:

Array of shape (6,) containing the shape function values at the given reference coordinates.

Return type:

numpy.ndarray

quadrature(nqp)

Return the quadrature weights and points for the specified number of quadrature points. Supported values for nqp are:

Parameters:

nqp (int) – Number of quadrature points. Supported values are 1 and 6.

Returns:

Tuple containing the quadrature weights and points.

Return type:

tuple(numpy.ndarray, numpy.ndarray)

quadrature_1pt()

1 point quadrature rule for a wedge element in 3D, which corresponds to evaluating the integrand at the centroid of the reference element.

\[\begin{split}w_0 = 1, \quad \boldsymbol \xi_0 = \begin{bmatrix} 1/3 \\ 1/3 \\ 0 \end{bmatrix}\end{split}\]

with \(w_0\) the quadrature weight and \(\boldsymbol \xi_0\) the quadrature point.

Returns:

Tuple containing the quadrature weights and points. The weight is returned as a 1D array of shape (1,), and the points are returned as a 2D array of shape (1, 3).

Return type:

tuple(numpy.ndarray, numpy.ndarray)

quadrature_6pt()

\(3 \times 2\) point quadrature rule for a wedge element in 3D, which corresponds to the tensor product of the 3 point quadrature rule for a triangle and the 2 Gauss point quadrature rule for a line.

\[\begin{split}w_{ij} = w_i^{\triangle} \cdot w_j^{\text{line}}, \quad \boldsymbol \xi_{ij} = \begin{bmatrix} \xi_i^{\triangle} \\ \eta_i^{\triangle} \\ \zeta_j^{\text{line}} \end{bmatrix}\end{split}\]
Returns:

Tuple containing the quadrature weights and points. The weights are returned as a 1D array of shape (6,), and the points are returned as a 2D array of shape (6, 3).

Return type:

tuple(numpy.ndarray, numpy.ndarray)