Home Installation Demos C++ docs Python docs
A finite element. More...
#include <finite-element.h>
Public Types | |
using | scalar_type = F |
Scalar type. | |
Public Member Functions | |
FiniteElement (element::family family, cell::type cell_type, polyset::type poly_type, int degree, const std::vector< std::size_t > &value_shape, mdspan_t< const F, 2 > wcoeffs, const std::array< std::vector< mdspan_t< const F, 2 >>, 4 > &x, const std::array< std::vector< mdspan_t< const F, 4 >>, 4 > &M, int interpolation_nderivs, maps::type map_type, sobolev::space sobolev_space, bool discontinuous, int embedded_subdegree, int embedded_superdegree, element::lagrange_variant lvariant, element::dpc_variant dvariant, std::vector< int > dof_ordering={}) | |
Construct a finite element. More... | |
FiniteElement (const FiniteElement &element)=default | |
Copy constructor. | |
FiniteElement (FiniteElement &&element)=default | |
Move constructor. | |
~FiniteElement ()=default | |
Destructor. | |
FiniteElement & | operator= (const FiniteElement &element)=default |
Assignment operator. | |
FiniteElement & | operator= (FiniteElement &&element)=default |
Move assignment operator. | |
bool | operator== (const FiniteElement &e) const |
Check if two elements are the same. More... | |
std::size_t | hash () const |
Get a unique hash of this element. | |
std::array< std::size_t, 4 > | tabulate_shape (std::size_t nd, std::size_t num_points) const |
Array shape for tabulate basis values and derivatives at set of points. More... | |
std::pair< std::vector< F >, std::array< std::size_t, 4 > > | tabulate (int nd, impl::mdspan_t< const F, 2 > x) const |
Compute basis values and derivatives at set of points. More... | |
std::pair< std::vector< F >, std::array< std::size_t, 4 > > | tabulate (int nd, std::span< const F > x, std::array< std::size_t, 2 > shape) const |
Compute basis values and derivatives at set of points. More... | |
void | tabulate (int nd, impl::mdspan_t< const F, 2 > x, mdspan_t< F, 4 > basis) const |
Compute basis values and derivatives at set of points. More... | |
void | tabulate (int nd, std::span< const F > x, std::array< std::size_t, 2 > xshape, std::span< F > basis) const |
Compute basis values and derivatives at set of points. More... | |
cell::type | cell_type () const |
Get the element cell type. More... | |
polyset::type | polyset_type () const |
Get the element polyset type. More... | |
int | degree () const |
Get the element polynomial degree. More... | |
int | embedded_superdegree () const |
Lowest degree n such that the highest degree polynomial in this element is contained in a Lagrange (or vector Lagrange) element of degree n . More... | |
int | embedded_subdegree () const |
Highest degree n such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this element. More... | |
const std::vector< std::size_t > & | value_shape () const |
Element value tensor shape. More... | |
int | dim () const |
Dimension of the finite element space. More... | |
element::family | family () const |
The finite element family. More... | |
element::lagrange_variant | lagrange_variant () const |
Lagrange variant of the element. More... | |
element::dpc_variant | dpc_variant () const |
DPC variant of the element. More... | |
maps::type | map_type () const |
Map type for the element. More... | |
sobolev::space | sobolev_space () const |
Underlying Sobolev space for this element. More... | |
bool | discontinuous () const |
Indicates whether this element is the discontinuous variant. More... | |
bool | dof_transformations_are_permutations () const |
Indicates if the degree-of-freedom transformations are all permutations. | |
bool | dof_transformations_are_identity () const |
Indicates is the dof transformations are all the identity. | |
std::pair< std::vector< F >, std::array< std::size_t, 3 > > | push_forward (impl::mdspan_t< const F, 3 > U, impl::mdspan_t< const F, 3 > J, std::span< const F > detJ, impl::mdspan_t< const F, 3 > K) const |
Map function values from the reference to a physical cell. More... | |
std::pair< std::vector< F >, std::array< std::size_t, 3 > > | pull_back (impl::mdspan_t< const F, 3 > u, impl::mdspan_t< const F, 3 > J, std::span< const F > detJ, impl::mdspan_t< const F, 3 > K) const |
Map function values from a physical cell to the reference. More... | |
template<typename O , typename P , typename Q , typename R > | |
std::function< void(O &, const P &, const Q &, F, const R &)> | map_fn () const |
Return a function that performs the appropriate push-forward/pull-back for the element type. More... | |
const std::vector< std::vector< std::vector< int > > > & | entity_dofs () const |
Get the dofs on each topological entity: (vertices, edges, faces, cell) in that order. More... | |
const std::vector< std::vector< std::vector< int > > > & | entity_closure_dofs () const |
Get the dofs on the closure of each topological entity: (vertices, edges, faces, cell) in that order. More... | |
std::pair< std::vector< F >, std::array< std::size_t, 3 > > | base_transformations () const |
Get the base transformations. More... | |
std::map< cell::type, std::pair< std::vector< F >, std::array< std::size_t, 3 > > > | entity_transformations () const |
Return the entity dof transformation matrices. More... | |
void | permute (std::span< std::int32_t > d, std::uint32_t cell_info) const |
Permute indices associated with degree-of-freedoms on the reference element ordering to the globally consistent physical element degree-of-freedom ordering. More... | |
void | permute_inv (std::span< std::int32_t > d, std::uint32_t cell_info) const |
Perform the inverse of the operation applied by permute(). More... | |
void | permute_subentity_closure (std::span< std::int32_t > d, std::uint32_t cell_info, cell::type entity_type, int entity_index) const |
Permute indices associated with degree-of-freedoms on the closure of a sub-entity of the reference element. More... | |
void | permute_subentity_closure_inv (std::span< std::int32_t > d, std::uint32_t cell_info, cell::type entity_type, int entity_index) const |
Perform the inverse of the operation applied by permute_subentity_closure(). More... | |
void | permute_subentity_closure (std::span< std::int32_t > d, std::uint32_t entity_info, cell::type entity_type) const |
Permute indices associated with degree-of-freedoms on the closure of a sub-entity of the reference element. More... | |
void | permute_subentity_closure_inv (std::span< std::int32_t > d, std::uint32_t entity_info, cell::type entity_type) const |
Perform the inverse of the operation applied by permute_subentity_closure(). More... | |
template<typename T > | |
void | T_apply (std::span< T > u, int n, std::uint32_t cell_info) const |
Transform basis functions from the reference element ordering and orientation to the globally consistent physical element ordering and orientation. More... | |
template<typename T > | |
void | Tt_apply (std::span< T > u, int n, std::uint32_t cell_info) const |
Apply the transpose of the operator applied by T_apply(). More... | |
template<typename T > | |
void | Tt_inv_apply (std::span< T > u, int n, std::uint32_t cell_info) const |
Apply the inverse transpose of the operator applied by T_apply(). More... | |
template<typename T > | |
void | Tinv_apply (std::span< T > u, int n, std::uint32_t cell_info) const |
Apply the inverse of the operator applied by T_apply(). More... | |
template<typename T > | |
void | Tt_apply_right (std::span< T > u, int n, std::uint32_t cell_info) const |
Right(post)-apply the transpose of the operator applied by T_apply(). More... | |
template<typename T > | |
void | T_apply_right (std::span< T > u, int n, std::uint32_t cell_info) const |
Right(post)-apply the operator applied by T_apply(). More... | |
template<typename T > | |
void | Tinv_apply_right (std::span< T > u, int n, std::uint32_t cell_info) const |
Right(post)-apply the inverse of the operator applied by T_apply(). More... | |
template<typename T > | |
void | Tt_inv_apply_right (std::span< T > u, int n, std::uint32_t cell_info) const |
Right(post)-apply the transpose inverse of the operator applied by T_apply(). More... | |
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & | points () const |
Return the interpolation points. More... | |
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & | interpolation_matrix () const |
Return a matrix of weights interpolation. More... | |
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & | dual_matrix () const |
Get the dual matrix. More... | |
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & | wcoeffs () const |
Get the coefficients that define the polynomial set in terms of the orthonormal polynomials. More... | |
const std::array< std::vector< std::pair< std::vector< F >, std::array< std::size_t, 2 > > >, 4 > & | x () const |
Get the interpolation points for each subentity. More... | |
const std::array< std::vector< std::pair< std::vector< F >, std::array< std::size_t, 4 > > >, 4 > & | M () const |
Get the interpolation matrices for each subentity. More... | |
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & | coefficient_matrix () const |
Get the matrix of coefficients. More... | |
bool | has_tensor_product_factorisation () const |
Indicates whether or not this element can be represented as a product of elements defined on lower-dimensional reference cells. More... | |
std::vector< std::vector< FiniteElement< F > > > | get_tensor_product_representation () const |
Get the tensor product representation of this element. More... | |
bool | interpolation_is_identity () const |
Indicates whether or not the interpolation matrix for this element is an identity matrix. More... | |
int | interpolation_nderivs () const |
The number of derivatives needed when interpolating. | |
const std::vector< int > & | dof_ordering () const |
Get dof layout. | |
Detailed Description
template<std::floating_point F>
class basix::FiniteElement< F >
A finite element.
The basis of a finite element is stored as a set of coefficients, which are applied to the underlying expansion set for that cell type, when tabulating.
Constructor & Destructor Documentation
◆ FiniteElement()
basix::FiniteElement< F >::FiniteElement | ( | element::family | family, |
cell::type | cell_type, | ||
polyset::type | poly_type, | ||
int | degree, | ||
const std::vector< std::size_t > & | value_shape, | ||
mdspan_t< const F, 2 > | wcoeffs, | ||
const std::array< std::vector< mdspan_t< const F, 2 >>, 4 > & | x, | ||
const std::array< std::vector< mdspan_t< const F, 4 >>, 4 > & | M, | ||
int | interpolation_nderivs, | ||
maps::type | map_type, | ||
sobolev::space | sobolev_space, | ||
bool | discontinuous, | ||
int | embedded_subdegree, | ||
int | embedded_superdegree, | ||
element::lagrange_variant | lvariant, | ||
element::dpc_variant | dvariant, | ||
std::vector< int > | dof_ordering = {} |
||
) |
Construct a finite element.
Initialising a finite element calculates the basis functions of the finite element, in terms of the polynomial basis.
The below explanation uses Einstein notation.
The basis functions \({\phi_i}\) of a finite element are represented as a linear combination of polynomials \(\{p_j\}\) in an underlying polynomial basis that span the space of all d-dimensional polynomials up to order \(k \ (P_k^d)\):
\[ \phi_i = c_{ij} p_j \]
In some cases, the basis functions \(\{\phi_i\}\) do not span the full space \(P_k\), in which case we denote space spanned by the basis functions by \(\{q_k\}\), which can be represented by:
\[ q_i = b_{ij} p_j. \]
This leads to
\[ \phi_i = c^{\prime}_{ij} q_j = c^{\prime}_{ij} b_{jk} p_k, \]
and in matrix form:
\[ \phi = C^{\prime} B p \]
If the basis functions span the full space, then \( B \) is simply the identity.
The basis functions \(\phi_i\) are defined by a dual set of functionals \(\{f_i\}\). The basis functions are the functions in span{ \(q_k\)} such that
\[ f_i(\phi_j) = \delta_{ij} \]
and inserting the expression for \(\phi_{j}\):
\[ f_i(c^{\prime}_{jk}b_{kl}p_{l}) = c^{\prime}_{jk} b_{kl} f_i \left( p_{l} \right) \]
Defining a matrix D given by applying the functionals to each polynomial \(p_j\):
\[ [D] = d_{ij},\mbox{ where } d_{ij} = f_i(p_j), \]
we have:
\[ C^{\prime} B D^{T} = I \]
and
\[ C^{\prime} = (B D^{T})^{-1}. \]
Recalling that \(C = C^{\prime} B\), where \(C\) is the matrix form of \(c_{ij}\),
\[ C = (B D^{T})^{-1} B \]
This function takes the matrices \(B\) (wcoeffs
) and \(D\) (M
) as inputs and will internally compute \(C\).
The matrix \(BD^{T}\) can be obtained from an element by using dual_matrix(). The matrix \(C\) can be obtained from an element by using coefficient_matrix().
Example: Order 1 Lagrange elements on a triangle
On a triangle, the scalar expansion basis is:
\[ p_0 = \sqrt{2}/2 \qquad p_1 = \sqrt{3}(2x + y - 1) \qquad p_2 = 3y - 1 \]
These span the space \(P_1\).
Lagrange order 1 elements span the space P_1, so in this example, B (span_coeffs) is the identity matrix:
\[ B = \begin{bmatrix} 1 & 0 & 0 \\ 0 & 1 & 0 \\ 0 & 0 & 1 \end{bmatrix} \]
The functionals defining the Lagrange order 1 space are point evaluations at the three vertices of the triangle. The matrix D (dual) given by applying these to p_0 to p_2 is:
\[ \mbox{dual} = \begin{bmatrix} \sqrt{2}/2 & -\sqrt{3} & -1 \\ \sqrt{2}/2 & \sqrt{3} & -1 \\ \sqrt{2}/2 & 0 & 2 \end{bmatrix} \]
For this example, this function outputs the matrix:
\[ C = \begin{bmatrix} \sqrt{2}/3 & -\sqrt{3}/6 & -1/6 \\ \sqrt{2}/3 & \sqrt{3}/6 & -1/6 \\ \sqrt{2}/3 & 0 & 1/3 \end{bmatrix} \]
The basis functions of the finite element can be obtained by applying the matrix C to the vector \([p_0, p_1, p_2]\), giving:
\[ \begin{bmatrix} 1 - x - y \\ x \\ y \end{bmatrix} \]
Example: Order 1 Raviart-Thomas on a triangle
On a triangle, the 2D vector expansion basis is:
\[ \begin{matrix} p_0 & = & (\sqrt{2}/2, 0) \\ p_1 & = & (\sqrt{3}(2x + y - 1), 0) \\ p_2 & = & (3y - 1, 0) \\ p_3 & = & (0, \sqrt{2}/2) \\ p_4 & = & (0, \sqrt{3}(2x + y - 1)) \\ p_5 & = & (0, 3y - 1) \end{matrix} \]
These span the space \( P_1^2 \).
Raviart-Thomas order 1 elements span a space smaller than \( P_1^2 \), so B (span_coeffs) is not the identity. It is given by:
\[ B = \begin{bmatrix} 1 & 0 & 0 & 0 & 0 & 0 \\ 0 & 0 & 0 & 1 & 0 & 0 \\ 1/12 & \sqrt{6}/48 & -\sqrt{2}/48 & 1/12 & 0 & \sqrt{2}/24 \end{bmatrix} \]
Applying the matrix B to the vector \([p_0, p_1, ..., p_5]\) gives the basis of the polynomial space for Raviart-Thomas:
\[ \begin{bmatrix} \sqrt{2}/2 & 0 \\ 0 & \sqrt{2}/2 \\ \sqrt{2}x/8 & \sqrt{2}y/8 \end{bmatrix} \]
The functionals defining the Raviart-Thomas order 1 space are integral of the normal components along each edge. The matrix D (dual) given by applying these to \(p_0\) to \(p_5\) is:
\[ D = \begin{bmatrix} -\sqrt{2}/2 & -\sqrt{3}/2 & -1/2 & -\sqrt{2}/2 & -\sqrt{3}/2 & -1/2 \\ -\sqrt{2}/2 & \sqrt{3}/2 & -1/2 & 0 & 0 & 0 \\ 0 & 0 & 0 & \sqrt{2}/2 & 0 & -1 \end{bmatrix} \]
In this example, this function outputs the matrix:
\[ C = \begin{bmatrix} -\sqrt{2}/2 & -\sqrt{3}/2 & -1/2 & -\sqrt{2}/2 & -\sqrt{3}/2 & -1/2 \\ -\sqrt{2}/2 & \sqrt{3}/2 & -1/2 & 0 & 0 & 0 \\ 0 & 0 & 0 & \sqrt{2}/2 & 0 & -1 \end{bmatrix} \]
The basis functions of the finite element can be obtained by applying the matrix C to the vector \([p_0, p_1, ..., p_5]\), giving:
\[ \begin{bmatrix} -x & -y \\ x - 1 & y \\ -x & 1 - y \end{bmatrix} \]
- Parameters
-
[in] family The element family [in] cell_type The cell type [in] poly_type The polyset type [in] degree The degree of the element [in] interpolation_nderivs The number of derivatives that need to be used during interpolation [in] value_shape The value shape of the element [in] wcoeffs Matrices for the kth value index containing the expansion coefficients defining a polynomial basis spanning the polynomial space for this element. Shape is (dim(finite element polyset), dim(Legendre polynomials)) [in] x Interpolation points. Indices are (tdim, entity index, point index, dim) [in] M The interpolation matrices. Indices are (tdim, entity index, dof, vs, point_index, derivative) [in] map_type The type of map to be used to map values from the reference to a cell [in] sobolev_space The underlying Sobolev space for the element [in] discontinuous Indicates whether or not this is the discontinuous version of the element [in] embedded_subdegree The highest degree n such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this element [in] embedded_superdegree The highest degree n such that at least one polynomial of degree n is included in this element's polymonial set [in] lvariant The Lagrange variant of the element [in] dvariant The DPC variant of the element [in] dof_ordering DOF reordering: a mapping from the reference order to a new permuted order
Member Function Documentation
◆ operator==()
bool FiniteElement::operator== | ( | const FiniteElement< F > & | e | ) | const |
Check if two elements are the same.
- Note
- This operator compares the element properties, e.g. family, degree, etc, and not computed numerical data
- Returns
- True if elements are the same
◆ tabulate_shape()
|
inline |
Array shape for tabulate basis values and derivatives at set of points.
- Parameters
-
[in] nd The order of derivatives, up to and including, to compute. Use 0 for the basis functions only. [in] num_points Number of points that basis will be computed at.
- Returns
- The shape of the array to will filled when passed to tabulate().
◆ tabulate() [1/4]
std::pair< std::vector< F >, std::array< std::size_t, 4 > > FiniteElement::tabulate | ( | int | nd, |
impl::mdspan_t< const F, 2 > | x | ||
) | const |
Compute basis values and derivatives at set of points.
- Note
- The version of tabulate() with the basis data as an out argument should be preferred for repeated call where performance is critical.
- Parameters
-
[in] nd The order of derivatives, up to and including, to compute. Use 0 for the basis functions only. [in] x The points at which to compute the basis functions. The shape of x is (number of points, geometric dimension).
- Returns
- The basis functions (and derivatives). The shape is (derivative, point, basis fn index, value index).
- The first index is the derivative, with higher derivatives are stored in triangular (2D) or tetrahedral (3D) ordering, ie for the (x,y) derivatives in 2D: (0,0), (1,0), (0,1), (2,0), (1,1), (0,2), (3,0)... The function basix::indexing::idx can be used to find the appropriate derivative.
- The second index is the point index
- The third index is the basis function index
- The fourth index is the basis function component. Its has size one for scalar basis functions.
◆ tabulate() [2/4]
std::pair< std::vector< F >, std::array< std::size_t, 4 > > FiniteElement::tabulate | ( | int | nd, |
std::span< const F > | x, | ||
std::array< std::size_t, 2 > | shape | ||
) | const |
Compute basis values and derivatives at set of points.
- Note
- The version of tabulate() with the basis data as an out argument should be preferred for repeated call where performance is critical
- Parameters
-
[in] nd The order of derivatives, up to and including, to compute. Use 0 for the basis functions only. [in] x The points at which to compute the basis functions (row-major storage). [in] shape The shape (number of points, geometric dimension)
of thex
array.
- Returns
- The basis functions (and derivatives). The shape is (derivative, point, basis fn index, value index).
- The first index is the derivative, with higher derivatives are stored in triangular (2D) or tetrahedral (3D) ordering, ie for the (x,y) derivatives in 2D: (0,0), (1,0), (0,1), (2,0), (1,1), (0,2), (3,0)... The function indexing::idx can be used to find the appropriate derivative.
- The second index is the point index
- The third index is the basis function index
- The fourth index is the basis function component. Its has size one for scalar basis functions.
◆ tabulate() [3/4]
void FiniteElement::tabulate | ( | int | nd, |
impl::mdspan_t< const F, 2 > | x, | ||
mdspan_t< F, 4 > | basis | ||
) | const |
Compute basis values and derivatives at set of points.
- Note
- This function is designed to be called at runtime, so its performance is critical.
- Parameters
-
[in] nd The order of derivatives, up to and including, to compute. Use 0 for the basis functions only. [in] x The points at which to compute the basis functions. The shape of x is (number of points, geometric dimension). [out] basis Memory location to fill. It must be allocated with shape (num_derivatives, num_points, num basis functions, value_size)
. The function tabulate_shape() can be used to get the required shape.- The first index is the derivative, with higher derivatives are stored in triangular (2D) or tetrahedral (3D) ordering, ie for the (x,y) derivatives in 2D: (0,0), (1,0), (0,1), (2,0), (1,1), (0,2), (3,0)... The function indexing::idx can be used to find the appropriate derivative.
- The second index is the point index
- The third index is the basis function index
- The fourth index is the basis function component. Its has size one for scalar basis functions.
- Todo:
- Remove all internal dynamic memory allocation, pass scratch space as required
◆ tabulate() [4/4]
void FiniteElement::tabulate | ( | int | nd, |
std::span< const F > | x, | ||
std::array< std::size_t, 2 > | xshape, | ||
std::span< F > | basis | ||
) | const |
Compute basis values and derivatives at set of points.
- Note
- This function is designed to be called at runtime, so its performance is critical.
- Parameters
-
[in] nd The order of derivatives, up to and including, to compute. Use 0 for the basis functions only. [in] x The points at which to compute the basis functions (row-major storage). The shape of x
is(number of points, geometric dimension)
.[in] xshape The shape (number of points, geometric dimension)
ofx
.[out] basis Memory location to fill. It must be allocated with shape (num_derivatives, num_points, num basis functions, value_size)
. The function tabulate_shape() can be used to get the required shape.- The first index is the derivative, with higher derivatives are stored in triangular (2D) or tetrahedral (3D) ordering, ie for the (x,y) derivatives in 2D: (0,0), (1,0), (0,1), (2,0), (1,1), (0,2), (3,0)... The function indexing::idx can be used to find the appropriate derivative.
- The second index is the point index
- The third index is the basis function index
- The fourth index is the basis function component. Its has size one for scalar basis functions.
◆ cell_type()
|
inline |
Get the element cell type.
- Returns
- The cell type
◆ polyset_type()
|
inline |
Get the element polyset type.
- Returns
- The polyset
◆ degree()
|
inline |
Get the element polynomial degree.
- Returns
- Polynomial degree
◆ embedded_superdegree()
|
inline |
Lowest degree n
such that the highest degree polynomial in this element is contained in a Lagrange (or vector Lagrange) element of degree n
.
- Returns
- Polynomial degree
◆ embedded_subdegree()
|
inline |
Highest degree n
such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this element.
- Returns
- Polynomial degree
◆ value_shape()
|
inline |
Element value tensor shape.
For example, returns {}
for scalars, {3}
for vectors in 3D, {2, 2}
for a rank-2 tensor in 2D.
- Returns
- Value shape
◆ dim()
|
inline |
Dimension of the finite element space.
The dimension is the number of degrees-of-freedom for the element.
- Returns
- Number of degrees of freedom
◆ family()
|
inline |
The finite element family.
- Returns
- The family
◆ lagrange_variant()
|
inline |
Lagrange variant of the element.
- Returns
- The Lagrange variant.
◆ dpc_variant()
|
inline |
DPC variant of the element.
- Returns
- The DPC variant.
◆ map_type()
|
inline |
Map type for the element.
- Returns
- The map type.
◆ sobolev_space()
|
inline |
Underlying Sobolev space for this element.
- Returns
- The Sobolev space.
◆ discontinuous()
|
inline |
Indicates whether this element is the discontinuous variant.
- Returns
- True if this element is a discontinuous version of the element.
◆ push_forward()
std::pair< std::vector< F >, std::array< std::size_t, 3 > > FiniteElement::push_forward | ( | impl::mdspan_t< const F, 3 > | U, |
impl::mdspan_t< const F, 3 > | J, | ||
std::span< const F > | detJ, | ||
impl::mdspan_t< const F, 3 > | K | ||
) | const |
Map function values from the reference to a physical cell.
This function can perform the mapping for multiple points, grouped by points that share a common Jacobian.
- Parameters
-
[in] U The function values on the reference. The indices are [Jacobian index, point index, components]
.[in] J The Jacobian of the mapping. The indices are [Jacobian index, J_i, J_j]
.[in] detJ The determinant of the Jacobian of the mapping. It has length J.shape(0)
[in] K The inverse of the Jacobian of the mapping. The indices are [Jacobian index, K_i, K_j]
.
- Returns
- The function values on the cell. The indices are [Jacobian index, point index, components].
◆ pull_back()
std::pair< std::vector< F >, std::array< std::size_t, 3 > > FiniteElement::pull_back | ( | impl::mdspan_t< const F, 3 > | u, |
impl::mdspan_t< const F, 3 > | J, | ||
std::span< const F > | detJ, | ||
impl::mdspan_t< const F, 3 > | K | ||
) | const |
Map function values from a physical cell to the reference.
- Parameters
-
[in] u The function values on the cell [in] J The Jacobian of the mapping [in] detJ The determinant of the Jacobian of the mapping [in] K The inverse of the Jacobian of the mapping
- Returns
- The function values on the reference. The indices are [Jacobian index, point index, components].
◆ map_fn()
|
inline |
Return a function that performs the appropriate push-forward/pull-back for the element type.
- Template Parameters
-
O The type that hold the (computed) mapped data (ndim==2) P The type that hold the data to be mapped (ndim==2) Q The type that holds the Jacobian (or inverse) matrix (ndim==2) R The type that holds the inverse of the Q
data (ndim==2)
- Returns
- A function that for a push-forward takes arguments
u
[out] The data on the physical cell after the push-forward flattened with row-major layout, shape=(num_points, value_size)U
[in] The data on the reference cell physical field to push forward, flattened with row-major layout, shape=(num_points, ref_value_size)J
[in] The Jacobian matrix of the map ,shape=(gdim, tdim)detJ
[in] det(J)K
[in] The inverse of the Jacobian matrix, shape=(tdim, gdim)
For a pull-back the arguments should be:
U
[out] The data on the reference cell after the pull-back, flattened with row-major layout, shape=(num_points, ref value_size)u
[in] The data on the physical cell that should be pulled back , flattened with row-major layout, shape=(num_points, value_size)K
[in] The inverse of the Jacobian matrix of the map ,shape=(tdim, gdim)detJ_inv
[in] 1/det(J)J
[in] The Jacobian matrix, shape=(gdim, tdim)
◆ entity_dofs()
|
inline |
Get the dofs on each topological entity: (vertices, edges, faces, cell) in that order.
For example, Lagrange degree 2 on a triangle has vertices: [[0], [1], [2]], edges: [[3], [4], [5]], cell: [[]]
- Returns
- Dofs associated with an entity of a given topological dimension. The shape is (tdim + 1, num_entities, num_dofs).
◆ entity_closure_dofs()
|
inline |
Get the dofs on the closure of each topological entity: (vertices, edges, faces, cell) in that order.
For example, Lagrange degree 2 on a triangle has vertices: [[0], [1], [2]], edges: [[1, 2, 3], [0, 2, 4], [0, 1, 5]], cell: [[0, 1, 2, 3, 4, 5]]
- Returns
- Dofs associated with the closure of an entity of a given topological dimension. The shape is
(tdim + 1, num_entities, num_dofs)
.
◆ base_transformations()
std::pair< std::vector< F >, std::array< std::size_t, 3 > > FiniteElement::base_transformations |
Get the base transformations.
The base transformations represent the effect of rotating or reflecting a subentity of the cell on the numbering and orientation of the DOFs. This returns a list of matrices with one matrix for each subentity permutation in the following order: Reversing edge 0, reversing edge 1, ... Rotate face 0, reflect face 0, rotate face 1, reflect face 1, ...
Example: Order 3 Lagrange on a triangle
This space has 10 dofs arranged like:
For this element, the base transformations are: [Matrix swapping 3 and 4, Matrix swapping 5 and 6, Matrix swapping 7 and 8] The first row shows the effect of reversing the diagonal edge. The second row shows the effect of reversing the vertical edge. The third row shows the effect of reversing the horizontal edge.
Example: Order 1 Raviart-Thomas on a triangle
This space has 3 dofs arranged like:
These DOFs are integrals of normal components over the edges: DOFs 0 and 2 are oriented inward, DOF 1 is oriented outwards. For this element, the base transformation matrices are:
The first matrix reverses DOF 0 (as this is on the first edge). The second matrix reverses DOF 1 (as this is on the second edge). The third matrix reverses DOF 2 (as this is on the third edge).
Example: DOFs on the face of Order 2 Nedelec first kind on a tetrahedron
On a face of this tetrahedron, this space has two face tangent DOFs:
For these DOFs, the subblocks of the base transformation matrices are:
- Returns
- The base transformations for this element. The shape is (ntranformations, ndofs, ndofs)
◆ entity_transformations()
|
inline |
Return the entity dof transformation matrices.
- Returns
- The entity transformations for the sub-entities of this element. The shape for each cell is (ntransformations, ndofs, ndofs)
◆ permute()
|
inline |
Permute indices associated with degree-of-freedoms on the reference element ordering to the globally consistent physical element degree-of-freedom ordering.
Given an array \(\tilde{d}\) that holds an integer associated with each degree-of-freedom and following the reference element degree-of-freedom ordering, this function computes
\[ d = P \tilde{d}, \]
where \(P\) is a permutation matrix and \(d\) holds the integers in \(\tilde{d}\) but permuted to follow the globally consistent physical element degree-of-freedom ordering. The permutation is computed in-place.
- Note
- This function is designed to be called at runtime, so its performance is critical.
- Parameters
-
[in,out] d Indices associated with each reference element degree-of-freedom (in). Indices associated with each physical element degree-of-freedom (out). cell_info Permutation info for the cell
◆ permute_inv()
|
inline |
Perform the inverse of the operation applied by permute().
Given an array \(d\) that holds an integer associated with each degree-of-freedom and following the globally consistent physical element degree-of-freedom ordering, this function computes
\[ \tilde{d} = P^{T} d, \]
where \(P^{T}\) is a permutation matrix and \(\tilde{d}\) holds the integers in \(d\) but permuted to follow the reference element degree-of-freedom ordering. The permutation is computed in-place.
- Parameters
-
[in,out] d Indices associated with each physical element degree-of-freedom [in]. Indices associated with each reference element degree-of-freedom [out]. cell_info Permutation info for the cell
◆ permute_subentity_closure() [1/2]
|
inline |
Permute indices associated with degree-of-freedoms on the closure of a sub-entity of the reference element.
This function performs a similar permutation to permute() but additionally permutes the positions of vertices and edges
- Note
- This function is designed to be called at runtime, so its performance is critical.
- Parameters
-
[in,out] d Indices associated with each reference element degree-of-freedom (in). Indices associated with each physical element degree-of-freedom (out). cell_info Permutation info for the cell entity_type The cell type of the sub-entity entity_index The index of the entity
◆ permute_subentity_closure_inv() [1/2]
|
inline |
Perform the inverse of the operation applied by permute_subentity_closure().
- Note
- This function is designed to be called at runtime, so its performance is critical.
- Parameters
-
[in,out] d Indices associated with each reference element degree-of-freedom (in). Indices associated with each physical element degree-of-freedom (out). cell_info Permutation info for the cell entity_type The cell type of the sub-entity entity_index The index of the entity
◆ permute_subentity_closure() [2/2]
|
inline |
Permute indices associated with degree-of-freedoms on the closure of a sub-entity of the reference element.
This function performs a similar permutation to permute() but additionally permutes the positions of vertices and edges
- Note
- This function is designed to be called at runtime, so its performance is critical.
- Parameters
-
[in,out] d Indices associated with each reference element degree-of-freedom (in). Indices associated with each physical element degree-of-freedom (out). entity_info Permutation info for the entity entity_type The cell type of the sub-entity
◆ permute_subentity_closure_inv() [2/2]
|
inline |
Perform the inverse of the operation applied by permute_subentity_closure().
- Note
- This function is designed to be called at runtime, so its performance is critical.
- Parameters
-
[in,out] d Indices associated with each reference element degree-of-freedom (in). Indices associated with each physical element degree-of-freedom (out). entity_info Permutation info for the entity entity_type The cell type of the sub-entity
◆ T_apply()
void basix::FiniteElement< F >::T_apply | ( | std::span< T > | u, |
int | n, | ||
std::uint32_t | cell_info | ||
) | const |
Transform basis functions from the reference element ordering and orientation to the globally consistent physical element ordering and orientation.
Consider that the value of a finite element function \(f_{h}\) at a point is given by
\[ f_{h} = \phi^{T} c, \]
where \(f_{h}\) has shape \(r \times 1\), \(\phi\) has shape \(d \times r\) and holds the finite element basis functions, and \(c\) has shape \(d \times 1\) and holds the degrees-of-freedom. The basis functions and degree-of-freedom are with respect to the physical element orientation. If the degrees-of-freedom on the physical element orientation are given by
\[ \phi = T \tilde{\phi}, \]
where \(T\) is a \(d \times d\) matrix, it follows from \(f_{h} = \phi^{T} c = \tilde{\phi}^{T} T^{T} c\) that
\[ \tilde{c} = T^{T} c. \]
This function applies \(T\) to data. The transformation is performed in-place. The operator \(T\) is orthogonal for many elements, but not all.
- Parameters
-
[in,out] u Data to transform. The shape is (m, n)
, wherem
is the number of dgerees-of-freedom and the storage is row-major.[in] n Number of columns in data
.cell_info Permutation info for the cell
◆ Tt_apply()
void basix::FiniteElement< F >::Tt_apply | ( | std::span< T > | u, |
int | n, | ||
std::uint32_t | cell_info | ||
) | const |
Apply the transpose of the operator applied by T_apply().
The transformation
\[ u \leftarrow T^{T} u \]
is performed in-place.
- Parameters
-
[in,out] u Data to transform. The shape is (m, n)
, wherem
is the number of dgerees-of-freedom an d the storage is row-major.[in] n Number of columns in data
.[in] cell_info Permutation info for the cell,
◆ Tt_inv_apply()
void basix::FiniteElement< F >::Tt_inv_apply | ( | std::span< T > | u, |
int | n, | ||
std::uint32_t | cell_info | ||
) | const |
Apply the inverse transpose of the operator applied by T_apply().
The transformation
\[ u \leftarrow T^{-T} u \]
is performed in-place.
- Parameters
-
[in,out] u Data to transform. The shape is (m, n)
, wherem
is the number of dgerees-of-freedom and the storage is row-major.[in] n Number of columns in data
.[in] cell_info Permutation info for the cell.
◆ Tinv_apply()
void basix::FiniteElement< F >::Tinv_apply | ( | std::span< T > | u, |
int | n, | ||
std::uint32_t | cell_info | ||
) | const |
Apply the inverse of the operator applied by T_apply().
The transformation
\[ u \leftarrow T^{-1} u \]
is performed in-place.
- Parameters
-
[in,out] u Data to transform. The shape is (m, n)
, wherem
is the number of dgerees-of-freedom and the storage is row-major.[in] n Number of columns in data
.[in] cell_info Permutation info for the cell.
◆ Tt_apply_right()
void basix::FiniteElement< F >::Tt_apply_right | ( | std::span< T > | u, |
int | n, | ||
std::uint32_t | cell_info | ||
) | const |
Right(post)-apply the transpose of the operator applied by T_apply().
Computes
\[ u^{T} \leftarrow u^{T} T^{T} \]
in-place.
- Parameters
-
[in,out] u Data to transform. The shape is (m, n)
, wherem
is the number of dgerees-of-freedom and the storage is row-major.[in] n Number of columns in data
.[in] cell_info Permutation info for the cell.
◆ T_apply_right()
void basix::FiniteElement< F >::T_apply_right | ( | std::span< T > | u, |
int | n, | ||
std::uint32_t | cell_info | ||
) | const |
Right(post)-apply the operator applied by T_apply().
Computes
\[ u^{T} \leftarrow u^{T} T \]
in-place.
- Parameters
-
[in,out] u Data to transform. The shape is (m, n)
, wherem
is the number of dgerees-of-freedom and the storage is row-major.[in] n Number of columns in data
.[in] cell_info Permutation info for the cell.
◆ Tinv_apply_right()
void basix::FiniteElement< F >::Tinv_apply_right | ( | std::span< T > | u, |
int | n, | ||
std::uint32_t | cell_info | ||
) | const |
Right(post)-apply the inverse of the operator applied by T_apply().
Computes
\[ u^{T} \leftarrow u^{T} T^{-1} \]
in-place.
- Parameters
-
[in,out] u Data to transform. The shape is (m, n)
, wherem
is the number of dgerees-of-freedom and the storage is row-major.[in] n Number of columns in data
.cell_info Permutation info for the cell.
◆ Tt_inv_apply_right()
void basix::FiniteElement< F >::Tt_inv_apply_right | ( | std::span< T > | u, |
int | n, | ||
std::uint32_t | cell_info | ||
) | const |
Right(post)-apply the transpose inverse of the operator applied by T_apply().
Computes
\[ u^{T} \leftarrow u^{T} T^{-T} \]
in-place.
- Parameters
-
[in,out] u Data to transform. The shape is (m, n)
, wherem
is the number of dgerees-of-freedom and the storage is row-major.[in] n Number of columns in data
.cell_info Permutation info for the cell.
◆ points()
|
inline |
Return the interpolation points.
The interpolation points are the coordinates on the reference element where a function need to be evaluated in order to interpolate it in the finite element space.
- Returns
- Array of coordinate with shape
(num_points, tdim)
◆ interpolation_matrix()
|
inline |
Return a matrix of weights interpolation.
To interpolate a function in this finite element, the functions should be evaluated at each point given by points(). These function values should then be multiplied by the weight matrix to give the coefficients of the interpolated function.
The shape of the returned matrix will be (dim, num_points * value_size)
, where dim
is the number of DOFs in the finite element, num_points
is the number of points returned by points(), and value_size
is the value size of the finite element.
For example, to interpolate into a Lagrange space, the following should be done:
To interpolate into a Raviart-Thomas space, the following should be done:
To interpolate into a Lagrange space with a block size, the following should be done:
- Returns
- The interpolation matrix. Shape is
(ndofs, number of interpolation points)
.
◆ dual_matrix()
|
inline |
Get the dual matrix.
This is the matrix \(BD^{T}\), as described in the documentation of the FiniteElement() constructor.
◆ wcoeffs()
|
inline |
Get the coefficients that define the polynomial set in terms of the orthonormal polynomials.
The polynomials spanned by each finite element in Basix are represented as a linear combination of the orthonormal polynomials of a given degree on the cell. Each row of this matrix defines a polynomial in the set spanned by the finite element.
For example, the orthonormal polynomials of degree <= 1 on a triangle are (where a, b, c, d are some constants):
- (sqrt(2), 0)
- (a*x - b, 0)
- (c*y - d, 0)
- (0, sqrt(2))
- (0, a*x - b)
- (0, c*y - d)
For a degree 1 Raviart-Thomas element, the first two rows of wcoeffs would be the following, as (1, 0) and (0, 1) are spanned by the element
- [1, 0, 0, 0, 0, 0]
- [0, 0, 0, 1, 0, 0]
The third row of wcoeffs in this example would give coefficients that represent (x, y) in terms of the orthonormal polynomials:
- [-b/(a*sqrt(2)), 1/a, 0, -d/(c*sqrt(2)), 0, 1/c]
These coefficients are only stored for custom elements. This function will throw an exception if called on a non-custom element.
- Returns
- Coefficient matrix. Shape is
(dim(finite element polyset), dim(Lagrange polynomials))
.
◆ x()
|
inline |
Get the interpolation points for each subentity.
The indices of this data are (tdim, entity index, point index, dim)
.
◆ M()
|
inline |
Get the interpolation matrices for each subentity.
The shape of this data is (tdim, entity index, dof, value size, point_index, derivative)
.
These matrices define how to evaluate the DOF functionals associated with each sub-entity of the cell. Given a function f, the functionals associated with the e
-th entity of dimension d
can be computed as follows:
For example, for a degree 1 Raviart-Thomas (RT) element on a triangle, the DOF functionals are integrals over the edges of the dot product of the function with the normal to the edge. In this case, x()
would contain quadrature points for each edge, and M()
would be a 1 by 2 by npoints
by 1 array for each edge. For each point, the [0, :, point, 0]
slice of this would be the quadrature weight multiplied by the normal. For all entities that are not edges, the entries in x()
and M()
for a degree 1 RT element would have size 0.
These matrices are only stored for custom elements. This function will throw an exception if called on a non-custom element
- Returns
- The interpolation matrices. The indices of this data are
(tdim, entity index, dof, vs, point_index, derivative)
.
◆ coefficient_matrix()
|
inline |
Get the matrix of coefficients.
- Returns
- The coefficient matrix. Shape is
(ndofs, ndofs)
.
◆ has_tensor_product_factorisation()
|
inline |
Indicates whether or not this element can be represented as a product of elements defined on lower-dimensional reference cells.
If the product exists, this element's basis functions can be computed as a tensor product of the basis elements of the elements in the product.
If such a factorisation exists, get_tensor_product_representation() can be used to get these elements.
◆ get_tensor_product_representation()
|
inline |
Get the tensor product representation of this element.
- Exceptions
-
std::runtime_error Thrown if no such factorisation exists.
The tensor product representation will be a vector of vectors of finite elements. Each tuple contains a vector of finite elements, and a vector of integers. The vector of finite elements gives the elements on an interval that appear in the tensor product representation. The vector of integers gives the permutation between the numbering of the tensor product DOFs and the number of the DOFs of this Basix element.
- Returns
- The tensor product representation
◆ interpolation_is_identity()
|
inline |
Indicates whether or not the interpolation matrix for this element is an identity matrix.
- Returns
- True if the interpolation matrix is the identity and false otherwise.
The documentation for this class was generated from the following files:
- /home/runner/work/basix/basix/cpp/basix/finite-element.h
- /home/runner/work/basix/basix/cpp/basix/finite-element.cpp