Basix 0.5.1

Home     Installation     Demos     C++ docs     Python docs

Public Member Functions | List of all members
basix::FiniteElement Class Reference

A finite element. More...

#include <finite-element.h>

Public Member Functions

 FiniteElement (element::family family, cell::type cell_type, int degree, const std::vector< std::size_t > &value_shape, const cmdspan2_t &wcoeffs, const std::array< std::vector< cmdspan2_t >, 4 > &x, const std::array< std::vector< cmdspan4_t >, 4 > &M, int interpolation_nderivs, maps::type map_type, bool discontinuous, int highest_complete_degree, int highest_degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int >>> tensor_factors={})
 Construct a finite element. More...
 
 FiniteElement (element::family family, cell::type cell_type, int degree, const std::vector< std::size_t > &value_shape, const cmdspan2_t &wcoeffs, const std::array< std::vector< cmdspan2_t >, 4 > &x, const std::array< std::vector< cmdspan4_t >, 4 > &M, int interpolation_nderivs, maps::type map_type, bool discontinuous, int highest_complete_degree, int highest_degree, element::lagrange_variant lvariant, std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int >>> tensor_factors={})
 Overloaded constructor.
 
 FiniteElement (element::family family, cell::type cell_type, int degree, const std::vector< std::size_t > &value_shape, const cmdspan2_t &wcoeffs, const std::array< std::vector< cmdspan2_t >, 4 > &x, const std::array< std::vector< cmdspan4_t >, 4 > &M, int interpolation_nderivs, maps::type map_type, bool discontinuous, int highest_complete_degree, int highest_degree, element::dpc_variant dvariant, std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int >>> tensor_factors={})
 Overloaded constructor.
 
 FiniteElement (element::family family, cell::type cell_type, int degree, const std::vector< std::size_t > &value_shape, const cmdspan2_t &wcoeffs, const std::array< std::vector< cmdspan2_t >, 4 > &x, const std::array< std::vector< cmdspan4_t >, 4 > &M, int interpolation_nderivs, maps::type map_type, bool disccontinuous, int highest_complete_degree, int highest_degree, std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int >>> tensor_factors={})
 Overloaded constructor.
 
 FiniteElement (const FiniteElement &element)=default
 Copy constructor.
 
 FiniteElement (FiniteElement &&element)=default
 Move constructor.
 
 ~FiniteElement ()=default
 Destructor.
 
FiniteElementoperator= (const FiniteElement &element)=default
 Assignment operator.
 
FiniteElementoperator= (FiniteElement &&element)=default
 Move assignment operator.
 
bool operator== (const FiniteElement &e) const
 
std::array< std::size_t, 4 > tabulate_shape (std::size_t nd, std::size_t num_points) const
 
std::pair< std::vector< double >, std::array< std::size_t, 4 > > tabulate (int nd, impl::cmdspan2_t x) const
 Compute basis values and derivatives at set of points. More...
 
std::pair< std::vector< double >, std::array< std::size_t, 4 > > tabulate (int nd, const std::span< const double > &x, std::array< std::size_t, 2 > shape) const
 
void tabulate (int nd, impl::cmdspan2_t x, impl::mdspan4_t basis) const
 
void tabulate (int nd, const std::span< const double > &x, std::array< std::size_t, 2 > xshape, const std::span< double > &basis) const
 
cell::type cell_type () const
 
int degree () const
 
int highest_degree () const
 
int highest_complete_degree () const
 
const std::vector< std::size_t > & value_shape () const
 
int dim () const
 
element::family family () const
 
element::lagrange_variant lagrange_variant () const
 
element::dpc_variant dpc_variant () const
 
maps::type map_type () const
 
bool discontinuous () const
 
bool dof_transformations_are_permutations () const
 
bool dof_transformations_are_identity () const
 
std::pair< std::vector< double >, std::array< std::size_t, 3 > > push_forward (impl::cmdspan3_t U, impl::cmdspan3_t J, std::span< const double > detJ, impl::cmdspan3_t K) const
 
std::pair< std::vector< double >, std::array< std::size_t, 3 > > pull_back (impl::cmdspan3_t u, impl::cmdspan3_t J, std::span< const double > detJ, impl::cmdspan3_t K) const
 
template<typename O , typename P , typename Q , typename R >
std::function< void(O &, const P &, const Q &, double, const R &)> map_fn () const
 
const std::vector< std::vector< std::vector< int > > > & entity_dofs () const
 
const std::vector< std::vector< std::vector< int > > > & entity_closure_dofs () const
 
std::pair< std::vector< double >, std::array< std::size_t, 3 > > base_transformations () const
 Get the base transformations. More...
 
std::map< cell::type, std::pair< std::vector< double >, std::array< std::size_t, 3 > > > entity_transformations () const
 
void permute_dofs (const std::span< std::int32_t > &dofs, std::uint32_t cell_info) const
 
void unpermute_dofs (const std::span< std::int32_t > &dofs, std::uint32_t cell_info) const
 
template<typename T >
void apply_dof_transformation (const std::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_transpose_dof_transformation (const std::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_inverse_transpose_dof_transformation (const std::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_inverse_dof_transformation (const std::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_dof_transformation_to_transpose (const std::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_transpose_dof_transformation_to_transpose (const std::span< T > &data, int block_size, std::uint32_t cell_info) const
 
template<typename T >
void apply_inverse_transpose_dof_transformation_to_transpose (const std::span< T > &data, int block_size, std::uint32_t cell_info) const
 Apply inverse transpose DOF transformations to some transposed data. More...
 
template<typename T >
void apply_inverse_dof_transformation_to_transpose (const std::span< T > &data, int block_size, std::uint32_t cell_info) const
 
const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & points () const
 
const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & interpolation_matrix () const
 Return a matrix of weights interpolation,. More...
 
const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & dual_matrix () const
 
const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & wcoeffs () const
 
const std::array< std::vector< std::pair< std::vector< double >, std::array< std::size_t, 2 > > >, 4 > & x () const
 
const std::array< std::vector< std::pair< std::vector< double >, std::array< std::size_t, 4 > > >, 4 > & M () const
 
const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & coefficient_matrix () const
 
bool has_tensor_product_factorisation () const
 
std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int > > > get_tensor_product_representation () const
 
bool interpolation_is_identity () const
 
int interpolation_nderivs () const
 The number of derivatives needed when interpolating.
 

Detailed Description

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()

FiniteElement::FiniteElement ( element::family  family,
cell::type  cell_type,
int  degree,
const std::vector< std::size_t > &  value_shape,
const cmdspan2_t &  wcoeffs,
const std::array< std::vector< cmdspan2_t >, 4 > &  x,
const std::array< std::vector< cmdspan4_t >, 4 > &  M,
int  interpolation_nderivs,
maps::type  map_type,
bool  discontinuous,
int  highest_complete_degree,
int  highest_degree,
element::lagrange_variant  lvariant,
element::dpc_variant  dvariant,
std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int >>>  tensor_factors = {} 
)

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 the function dual_matrix(). The matrix \(C\) can be obtained from an element by using the function 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]familyThe element family
[in]cell_typeThe cell type
[in]degreeThe degree of the element
[in]interpolation_nderivsThe number of derivatives that need to be used during interpolation
[in]value_shapeThe value shape of the element
[in]wcoeffsMatrices 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]xInterpolation points. Indices are (tdim, entity index, point index, dim)
[in]MThe interpolation matrices. Indices are (tdim, entity index, dof, vs, point_index, derivative)
[in]map_typeThe type of map to be used to map values from the reference to a cell
[in]discontinuousIndicates whether or not this is the discontinuous version of the element
[in]highest_complete_degreeThe highest degree n such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this element
[in]highest_degreeThe highest degree n such that at least one polynomial of degree n is included in this element's polymonial set
[in]lvariantThe Lagrange variant of the element
[in]dvariantThe DPC variant of the element
[in]tensor_factorsThe factors in the tensor product representation of this element

Member Function Documentation

◆ apply_dof_transformation()

template<typename T >
void basix::FiniteElement::apply_dof_transformation ( const std::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply DOF transformations to some data

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_dof_transformation_to_transpose()

template<typename T >
void basix::FiniteElement::apply_dof_transformation_to_transpose ( const std::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply DOF transformations to some transposed data

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_inverse_dof_transformation()

template<typename T >
void basix::FiniteElement::apply_inverse_dof_transformation ( const std::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply inverse DOF transformations to some data

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_inverse_dof_transformation_to_transpose()

template<typename T >
void basix::FiniteElement::apply_inverse_dof_transformation_to_transpose ( const std::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply inverse DOF transformations to some transposed data

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_inverse_transpose_dof_transformation()

template<typename T >
void basix::FiniteElement::apply_inverse_transpose_dof_transformation ( const std::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply inverse transpose DOF transformations to some data

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_inverse_transpose_dof_transformation_to_transpose()

template<typename T >
void basix::FiniteElement::apply_inverse_transpose_dof_transformation_to_transpose ( const std::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply inverse transpose DOF transformations to some transposed data.

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_transpose_dof_transformation()

template<typename T >
void basix::FiniteElement::apply_transpose_dof_transformation ( const std::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply transpose DOF transformations to some data

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ apply_transpose_dof_transformation_to_transpose()

template<typename T >
void basix::FiniteElement::apply_transpose_dof_transformation_to_transpose ( const std::span< T > &  data,
int  block_size,
std::uint32_t  cell_info 
) const

Apply transpose DOF transformations to some transposed data

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dataThe data
block_sizeThe number of data points per DOF
cell_infoThe permutation info for the cell

◆ base_transformations()

std::pair< std::vector< double >, std::array< std::size_t, 3 > > FiniteElement::base_transformations ( ) const

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:

2
|\
6 4
| \
5 9 3
| \
0-7-8-1

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:

|\
| \
| \
<-1 0
| / \
| L ^ \
| | \
---2---

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:

0: [[-1, 0, 0],
[ 0, 1, 0],
[ 0, 0, 1]]
1: [[1, 0, 0],
[0, -1, 0],
[0, 0, 1]]
2: [[1, 0, 0],
[0, 1, 0],
[0, 0, -1]]

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:

|\ |\
| \ | \
| \ | ^\
| \ | | \
| 0->\ | 1 \
| \ | \
------ ------

For these DOFs, the subblocks of the base transformation matrices are:

rotation: [[-1, 1],
[ 1, 0]]
reflection: [[0, 1],
[1, 0]]
Returns
The base transformations for this element. The shape is (ntranformations, ndofs, ndofs)

◆ cell_type()

cell::type FiniteElement::cell_type ( ) const

Get the element cell type

Returns
The cell type

◆ coefficient_matrix()

const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & FiniteElement::coefficient_matrix ( ) const

Get the matrix of coefficients.

This is the matrix \(C\), as described in the documentation of the FiniteElement() constructor.

Returns
The coefficient matrix. Shape is (ndofs, ndofs)

◆ degree()

int FiniteElement::degree ( ) const

Get the element polynomial degree

Returns
Polynomial degree

◆ dim()

int FiniteElement::dim ( ) const

Dimension of the finite element space (number of degrees-of-freedom for the element)

Returns
Number of degrees of freedom

◆ discontinuous()

bool FiniteElement::discontinuous ( ) const

Indicates whether this element is the discontinuous variant

Returns
True if this element is a discontinuous version of the element

◆ dof_transformations_are_identity()

bool FiniteElement::dof_transformations_are_identity ( ) const

Indicates whether the dof transformations are all the identity

Returns
True or False

◆ dof_transformations_are_permutations()

bool FiniteElement::dof_transformations_are_permutations ( ) const

Indicates whether the dof transformations are all permutations

Returns
True or False

◆ dpc_variant()

element::dpc_variant FiniteElement::dpc_variant ( ) const

Get the DPC variant of the element.

Returns
The DPC variant

◆ dual_matrix()

const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & FiniteElement::dual_matrix ( ) const

Get the dual matrix.

This is the matrix \(BD^{T}\), as described in the documentation of the FiniteElement() constructor.

Returns
The dual matrix. Shape is (ndofs, ndofs)

◆ entity_closure_dofs()

const std::vector< std::vector< std::vector< int > > > & FiniteElement::entity_closure_dofs ( ) const

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).

◆ entity_dofs()

const std::vector< std::vector< std::vector< int > > > & FiniteElement::entity_dofs ( ) const

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_transformations()

std::map< cell::type, std::pair< std::vector< double >, std::array< std::size_t, 3 > > > FiniteElement::entity_transformations ( ) const

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)

◆ family()

element::family FiniteElement::family ( ) const

Get the finite element family

Returns
The family

◆ get_tensor_product_representation()

std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int > > > FiniteElement::get_tensor_product_representation ( ) const

Get the tensor product representation of this element, or throw an error if no such factorisation exists.

The tensor product representation will be a vector of tuples. 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

◆ has_tensor_product_factorisation()

bool FiniteElement::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. 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.

◆ highest_complete_degree()

int FiniteElement::highest_complete_degree ( ) const

Highest degree n such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this element.

Returns
Polynomial degree

◆ highest_degree()

int FiniteElement::highest_degree ( ) 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.

Returns
Polynomial degree

◆ interpolation_is_identity()

bool FiniteElement::interpolation_is_identity ( ) const

Indicates whether or not the interpolation matrix for this element is an identity matrix

◆ interpolation_matrix()

const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & FiniteElement::interpolation_matrix ( ) const

Return a matrix of weights interpolation,.

To interpolate a function in this finite element, the functions should be evaluated at each point given by FiniteElement::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:

i_m = element.interpolation_matrix()
pts = element.points()
values = vector(pts.shape(0))
FOR i, p IN ENUMERATE(pts):
values[i] = f.evaluate_at(p)
coefficients = i_m * values

To interpolate into a Raviart-Thomas space, the following should be done:

i_m = element.interpolation_matrix()
pts = element.points()
vs = prod(element.value_shape())
values = VECTOR(pts.shape(0) * vs)
FOR i, p IN ENUMERATE(pts):
values[i::pts.shape(0)] = f.evaluate_at(p)
coefficients = i_m * values

To interpolate into a Lagrange space with a block size, the following should be done:

i_m = element.interpolation_matrix()
pts = element.points()
coefficients = VECTOR(element.dim() * block_size)
FOR b IN RANGE(block_size):
values = vector(pts.shape(0))
FOR i, p IN ENUMERATE(pts):
values[i] = f.evaluate_at(p)[b]
coefficients[::block_size] = i_m * values
Returns
The interpolation matrix. Shape is (ndofs, number of interpolation points)

◆ lagrange_variant()

element::lagrange_variant FiniteElement::lagrange_variant ( ) const

Get the Lagrange variant of the element.

Returns
The Lagrange variant

◆ M()

const std::array< std::vector< std::pair< std::vector< double >, std::array< std::size_t, 4 > > >, 4 > & FiniteElement::M ( ) const

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:

matrix = element.M()[d][e]
pts = element.x()[d][e]
nderivs = element
values = f.eval_derivs(nderivs, pts)
result = ZEROS(matrix.shape(0))
FOR i IN RANGE(matrix.shape(0)):
FOR j IN RANGE(matrix.shape(1)):
FOR k IN RANGE(matrix.shape(2)):
FOR l IN RANGE(matrix.shape(3)):
result[i] += matrix[i, j, k, l] * values[l][k][j]

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 by 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)

◆ map_fn()

template<typename O , typename P , typename Q , typename R >
std::function<void(O&, const P&, const Q&, double, const R&)> basix::FiniteElement::map_fn ( ) const
inline

Return a function that performs the appropriate push-forward/pull-back for the element type

Template Parameters
OThe type that hold the (computed) mapped data (ndim==2)
PThe type that hold the data to be mapped (ndim==2)
QThe type that holds the Jacobian (or inverse) matrix (ndim==2)
RThe 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)

◆ map_type()

maps::type FiniteElement::map_type ( ) const

Get the map type for this element

Returns
The map type

◆ operator==()

bool FiniteElement::operator== ( const FiniteElement 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

◆ permute_dofs()

void FiniteElement::permute_dofs ( const std::span< std::int32_t > &  dofs,
std::uint32_t  cell_info 
) const

Permute the dof numbering on a cell

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dofsThe dof numbering for the cell
cell_infoThe permutation info for the cell

◆ points()

const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & FiniteElement::points ( ) const

Return the interpolation points, i.e. 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)

◆ pull_back()

std::pair< std::vector< double >, std::array< std::size_t, 3 > > FiniteElement::pull_back ( impl::cmdspan3_t  u,
impl::cmdspan3_t  J,
std::span< const double >  detJ,
impl::cmdspan3_t  K 
) const

Map function values from a physical cell to the reference

Parameters
[in]uThe function values on the cell
[in]JThe Jacobian of the mapping
[in]detJThe determinant of the Jacobian of the mapping
[in]KThe inverse of the Jacobian of the mapping
Returns
The function values on the reference. The indices are [Jacobian index, point index, components].

◆ push_forward()

std::pair< std::vector< double >, std::array< std::size_t, 3 > > FiniteElement::push_forward ( impl::cmdspan3_t  U,
impl::cmdspan3_t  J,
std::span< const double >  detJ,
impl::cmdspan3_t  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]UThe function values on the reference. The indices are [Jacobian index, point index, components].
[in]JThe Jacobian of the mapping. The indices are [Jacobian index, J_i, J_j].
[in]detJThe determinant of the Jacobian of the mapping. It has length J.shape(0)
[in]KThe 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].

◆ tabulate() [1/4]

std::pair< std::vector< double >, std::array< std::size_t, 4 > > FiniteElement::tabulate ( int  nd,
const std::span< const double > &  x,
std::array< std::size_t, 2 >  shape 
) const

Compute basis values and derivatives at set of points.

Note
The version of FiniteElement::tabulate with the basis data as an out argument should be preferred for repeated call where performance is critical
Parameters
[in]ndThe order of derivatives, up to and including, to compute. Use 0 for the basis functions only.
[in]xThe points at which to compute the basis functions (row-major storage).
[in]shapeThe shape (number of points, geometric dimension) of the x 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 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]

void FiniteElement::tabulate ( int  nd,
const std::span< const double > &  x,
std::array< std::size_t, 2 >  xshape,
const std::span< double > &  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]ndThe order of derivatives, up to and including, to compute. Use 0 for the basis functions only.
[in]xThe points at which to compute the basis functions (row-major storage). The shape of x is (number of points, geometric dimension).
[in]xshapeThe shape (number of points, geometric dimension) of x.
[out]basisMemory location to fill. It must be allocated with shape (num_derivatives, num_points, num basis functions, value_size). The function FiniteElement::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 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() [3/4]

std::pair< std::vector< double >, std::array< std::size_t, 4 > > FiniteElement::tabulate ( int  nd,
impl::cmdspan2_t  x 
) const

Compute basis values and derivatives at set of points.

Note
The version of FiniteElement::tabulate with the basis data as an out argument should be preferred for repeated call where performance is critical
Parameters
[in]ndThe order of derivatives, up to and including, to compute. Use 0 for the basis functions only.
[in]xThe 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() [4/4]

void FiniteElement::tabulate ( int  nd,
impl::cmdspan2_t  x,
impl::mdspan4_t  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]ndThe order of derivatives, up to and including, to compute. Use 0 for the basis functions only.
[in]xThe points at which to compute the basis functions. The shape of x is (number of points, geometric dimension).
[out]basisMemory location to fill. It must be allocated with shape (num_derivatives, num_points, num basis functions, value_size). The function FiniteElement::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 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.
Todo:
Remove all internal dynamic memory allocation, pass scratch space as required

◆ tabulate_shape()

std::array< std::size_t, 4 > FiniteElement::tabulate_shape ( std::size_t  nd,
std::size_t  num_points 
) const

Array shape for tabulate basis values and derivatives at set of points.

Parameters
[in]ndThe order of derivatives, up to and including, to compute. Use 0 for the basis functions only.
[in]num_pointsNumber of points that basis will be computed at
Returns
The shape of the array to will filled when passed to FiniteElement::tabulate

◆ unpermute_dofs()

void FiniteElement::unpermute_dofs ( const std::span< std::int32_t > &  dofs,
std::uint32_t  cell_info 
) const

Unpermute the dof numbering on a cell

Note
This function is designed to be called at runtime, so its performance is critical.
Parameters
[in,out]dofsThe dof numbering for the cell
cell_infoThe permutation info for the cell

◆ value_shape()

const std::vector< std::size_t > & FiniteElement::value_shape ( ) const

The element value tensor shape, e.g. returning {} for scalars, {3} for vectors in 3D, {2, 2} for a rank-2 tensor in 2D.

Returns
Value shape

◆ wcoeffs()

const std::pair< std::vector< double >, std::array< std::size_t, 2 > > & FiniteElement::wcoeffs ( ) const

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()

const std::array< std::vector< std::pair< std::vector< double >, std::array< std::size_t, 2 > > >, 4 > & FiniteElement::x ( ) const

Get the interpolation points for each subentity.

The indices of this data are (tdim, entity index, point index, dim).


The documentation for this class was generated from the following files: