elements

A mesh is composed of non-overlapping elements, where the elements sub-library contains the mathematical functions to support a very broad range of element types including:

  • linear, quadratic, and cubic serendipity elements in 2D and 3D;

  • arbitrary-order tensor product elements in 2D and 3D;

  • arbitrary-order spectral elements; and

  • a linear 4D element.

The elements sub-library has functions to calculate quantities that are commonly used in finite element methods (both continuous and discontinuous) such as a basis function, gradient of a basis function, the Jacobian matrix, the inverse Jacobian matrix, the determinant of the Jacobian matrix, and a physical position inside the element, to name a few examples. The elements sub-library also supports both Gauss-Legendre and Gauss-Lobatto quadrature rules up to 8 quadrature points in each coordinate direction.

elements

A mesh is a collection of non-overlapping elements, and the mesh can be unstructured. Each reference element is defined in a reference coordinate system with a spatial mapping (e.g., an interpolation polynomial) to the physical coordinates. An element in the physical coordinates can have edges that are defined by straight lines (i.e., linear) or have edges defined by a high-order polynomial. The elements sub-library contains the mathematical functions for a suite of element types. The geometry sub-library combines the mesh data structures and mesh index spaces in SWAGE with the reference element defined in elements.

geometry definition

The position inside an element is given by a polynomial. elements contains a suite of spatial polynomials and also a space-time polynomial. The Jacobi matrix is equal to the gradient of this spatial polynomial. The volume of an element can be calculated using the determinate of the Jacobi matrix. The inverse of the Jacobi matrix is used to transform a Nabla operator from the physical coordinates to the reference coordinates.

index conventions

The reference element contains nodes that are collocated with the Lobatto quadrature points. The local indexing of the quadrature points and nodes within an element follows an (i,j,k) index convention. An index map is supplied to convert the Serendipity local index convention to the local (i,j,k) index convention.

namespace elements

Functions

void lobatto_nodes_1D(CArray<real_t> &lob_nodes_1D, const int &num)

Elements contains many of the basis functions required to implement a wide range of numerical methods for computational physics. Currently 2D and 3D serendipity basis set is included, as well as the arbitrary order 2D and 3D tensor product hexahedral elements. There is also a 4D tesseract element that can be used for space-time methods. Make sure include the elements.h file in your code to access the elememnts library. lobatto_nodes_1D creates nodal positions defined on [-1,1] using the Gauss-Lobatto quadrature points. The CArray lob_nodes_1D is passed in by reference and modified in place. The integer num is the numer of points being defined in 1D.

void lobatto_weights_1D(CArray<real_t> &lob_weights_1D, const int &num)

lobatto_weights_1D creates quadrature weights corresponding the nodal positions defined on [-1,1] using the Gauss-Lobatto quadrature points. The CArray lob_weights_1D is passed in my reference and modified in place. The integer num is the numer of points being defined in 1D.

void legendre_nodes_1D(CArray<real_t> &leg_nodes_1D, const int &num)
void legendre_weights_1D(CArray<real_t> &leg_weights_1D, const int &num)
void length_weights(CArray<real_t> &len_weights_1D, CArray<real_t> &lob_weights_1D, CArray<real_t> &lob_nodes_1D, const int &p_order)

length_weights partitions quadrature weights corresponding the 1D nodal positions defined on [-1,1] using a distance weighted partition. The CArray len_weights_1D is passed in my reference and modified in place. The CArray lob_weights_1D and lob_nodes_1D are passed in and the weights are partitioned. The integer p_order is the numer of kinematic polynomial order of the desired element. Note: This function currently only supports up to third order elements.

void sub_weights(CArray<real_t> &sub_weights_1D, CArray<real_t> &lob_weights_1D, CArray<real_t> &lob_nodes_1D, const int &p_order)

sub_weights partitions quadrature weights to the corners corresponding the 1D nodal positions defined on [-1,1] using a distance weighted partition.

void mat_inverse(ViewCArray<real_t> &mat_inv, ViewCArray<real_t> &matrix)

mat_inverse is a light weight function for inverting a 3X3 matrix.

void set_nodes_wgts(CArray<real_t> &lob_nodes_1D, CArray<real_t> &lob_weights_1D, CArray<real_t> &len_weights_1D, CArray<real_t> &sub_weights_1D, const int p_order)

set_nodes_wgts takes in the CArrays used to hold the quadarture nodes, weights, and partitioned weights and initializes everyting based of the desired polynomial order

void set_unit_normals(CArray<real_t> &ref_corner_unit_normals)

set_unit_normals takes in the CArray defined to hold the unit normals of the corners in a cell in reference space and initializes them.

void line_gauss_info(real_t &x, real_t &w, int &m, int &p)

line_gauss_info is used to create the nodes and weights for Gauss-Legendre quadrature in 1D. This function is is used to create quadrature points in 2D, 3D, and 4D.

void line_lobatto_info(real_t &x, real_t &w, int &m, int &p)
void gauss_2d(ViewCArray<real_t> &these_g_pts, ViewCArray<real_t> &these_weights, ViewCArray<real_t> &tot_g_weight, int &quad_order)
void gauss_3d(CArray<real_t> &these_g_pts, CArray<real_t> &these_weights, CArray<real_t> &tot_g_weight, int &quad_order)
void gauss_4d(ViewCArray<real_t> &these_g_pts, ViewCArray<real_t> &these_weights, int &quad_order, const int &dim)
void lobatto_2d(ViewCArray<real_t> &these_L_pts, ViewCArray<real_t> &these_weights, int &quad_order)
void lobatto_3d(ViewCArray<real_t> &these_L_pts, ViewCArray<real_t> &these_weights, int &quad_order)
void lobatto_4d(ViewCArray<real_t> &these_L_pts, ViewCArray<real_t> &these_weights, int &quad_order, const int &dim)
void jacobian_2d(ViewCArray<real_t> &J_matrix, real_t &det_J, const ViewCArray<real_t> &vertices, const ViewCArray<real_t> &this_partial, const int &num_verts)
void jacobian_3d(ViewCArray<real_t> &J_matrix, real_t &det_J, const ViewCArray<real_t> &vertices, const ViewCArray<real_t> &this_partial, const int &num_verts)
void jacobian_4d(ViewCArray<real_t> &J_matrix, real_t &det_J, const ViewCArray<real_t> &vertices, const ViewCArray<real_t> &this_partial, const int &num_verts, const int &dim)
void jacobian_inverse_2d(ViewCArray<real_t> &J_inverse, const ViewCArray<real_t> &jacobian)
void jacobian_inverse_3d(ViewCArray<real_t> &J_inverse, const ViewCArray<real_t> &jacobian)
void jacobian_inverse_4d(ViewCArray<real_t> &J_inverse, const ViewCArray<real_t> &jacobian, const real_t &det_J)
void chebyshev_nodes_1D(ViewCArray<real_t> &cheb_nodes_1D, const int &order)
void gauss_3d(ViewCArray<real_t> &these_g_pts, ViewCArray<real_t> &these_weights, ViewCArray<real_t> &tot_g_weight, int &quad_order)
class Element2D

Subclassed by elements::Quad12, elements::Quad4, elements::Quad8

class Element3D

Subclassed by elements::Hex20, elements::Hex32, elements::Hex8

class Element4D

Subclassed by elements::Tess16

class Quad4 : public elements::Element2D
class Quad8 : public elements::Element2D
class Quad12 : public elements::Element2D
class QuadN
class Hex8 : public elements::Element3D
class Hex20 : public elements::Element3D
class Hex32 : public elements::Element3D
class HexN
class Tess16 : public elements::Element4D
struct elem_type_t
class element_selector
class ref_element
namespace elem_types

Enums

enum elem_type

Values:

enumerator Quad4
enumerator Quad8
enumerator Quad12
enumerator QuadN
enumerator Hex8
enumerator Hex20
enumerator Hex32
enumerator HexN