DOLFINx  0.4.1
DOLFINx C++ interface
Public Member Functions | List of all members
FiniteElement Class Reference

Finite Element, containing the dof layout on a reference element, and various methods for evaluating and transforming the basis. More...

#include <FiniteElement.h>

Public Member Functions

 FiniteElement (const ufcx_finite_element &e)
 Create finite element from UFC finite element. More...
 
 FiniteElement (const basix::FiniteElement &element, int bs)
 Create finite element from a Basix finite element. More...
 
 FiniteElement (const FiniteElement &element)=delete
 Copy constructor.
 
 FiniteElement (FiniteElement &&element)=default
 Move constructor.
 
virtual ~FiniteElement ()=default
 Destructor.
 
FiniteElementoperator= (const FiniteElement &element)=delete
 Copy assignment.
 
FiniteElementoperator= (FiniteElement &&element)=default
 Move assignment.
 
bool operator== (const FiniteElement &e) const
 Check if two elements are equivalent. More...
 
bool operator!= (const FiniteElement &e) const
 Check if two elements are not equivalent. More...
 
std::string signature () const noexcept
 String identifying the finite element. More...
 
mesh::CellType cell_shape () const noexcept
 Cell shape. More...
 
int space_dimension () const noexcept
 Dimension of the finite element function space (the number of degrees-of-freedom for the element) More...
 
int block_size () const noexcept
 Block size of the finite element function space. For VectorElements and TensorElements, this is the number of DOFs colocated at each DOF point. For other elements, this is always 1. More...
 
int value_size () const
 The value size, e.g. 1 for a scalar function, 2 for a 2D vector, 9 for a second-order tensor in 3D. More...
 
int reference_value_size () const
 The value size, e.g. 1 for a scalar function, 2 for a 2D vector, 9 for a second-order tensor in 3D, for the reference element. More...
 
xtl::span< const int > value_shape () const noexcept
 Shape of the value space. The rank is the size of the value_shape.
 
std::string family () const noexcept
 The finite element family. More...
 
void tabulate (xt::xtensor< double, 4 > &values, const xt::xtensor< double, 2 > &X, int order) const
 Evaluate all derivatives of the basis functions up to given order at given points in reference cell. More...
 
template<typename O , typename P , typename Q , typename R >
std::function< void(O &, const P &, const Q &, double, const R &)> map_fn () const
 Return a function that performs the appropriate push-forward (pull-back) for the element type. More...
 
int num_sub_elements () const noexcept
 Get the number of sub elements (for a mixed or blocked element) More...
 
bool is_mixed () const noexcept
 Check if element is a mixed element, i.e. composed of two or more elements of different types. A block element, e.g. a Lagrange element with block size > 1 is not considered mixed. More...
 
const std::vector< std::shared_ptr< const FiniteElement > > & sub_elements () const noexcept
 Subelements (if any)
 
std::shared_ptr< const FiniteElementextract_sub_element (const std::vector< int > &component) const
 Extract sub finite element for component.
 
basix::maps::type map_type () const
 Get the map type used by the element.
 
bool interpolation_ident () const noexcept
 Check if interpolation into the finite element space is an identity operation given the evaluation on an expression at specific points, i.e. the degree-of-freedom are equal to point evaluations. The function will return true for Lagrange elements. More...
 
bool map_ident () const noexcept
 Check if the push forward/pull back map from the values on reference to the values on a physical cell for this element is the identity map. More...
 
const xt::xtensor< double, 2 > & interpolation_points () const
 Points on the reference cell at which an expression need to be evaluated in order to interpolate the expression in the finite element space. For Lagrange elements the points will just be the nodal positions. For other elements the points will typically be the quadrature points used to evaluate moment degrees of freedom. More...
 
const xt::xtensor< double, 2 > & interpolation_operator () const
 Interpolation operator (matrix) Pi that maps a function evaluated at the points provided by FiniteElement::interpolation_points to the element degrees of freedom, i.e. dofs = Pi f_x. See the Basix documentation for basix::FiniteElement::interpolation_matrix for how the data in f_x should be ordered. More...
 
xt::xtensor< double, 2 > create_interpolation_operator (const FiniteElement &from) const
 Create a matrix that maps degrees of freedom from one element to this element (interpolation). More...
 
bool needs_dof_transformations () const noexcept
 Check if DOF transformations are needed for this element. More...
 
bool needs_dof_permutations () const noexcept
 Check if DOF permutations are needed for this element. More...
 
template<typename T >
std::function< void(const xtl::span< T > &, const xtl::span< const std::uint32_t > &, std::int32_t, int)> get_dof_transformation_function (bool inverse=false, bool transpose=false, bool scalar_element=false) const
 Return a function that applies DOF transformation to some data. More...
 
template<typename T >
std::function< void(const xtl::span< T > &, const xtl::span< const std::uint32_t > &, std::int32_t, int)> get_dof_transformation_to_transpose_function (bool inverse=false, bool transpose=false, bool scalar_element=false) const
 Return a function that applies DOF transformation to some transposed data. More...
 
template<typename T >
void apply_dof_transformation (const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
 Apply DOF transformation to some data. More...
 
template<typename T >
void apply_inverse_transpose_dof_transformation (const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
 Apply inverse transpose transformation to some data. For VectorElements, this applies the transformations for the scalar subelement. More...
 
template<typename T >
void apply_transpose_dof_transformation (const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
 Apply transpose transformation to some data. For VectorElements, this applies the transformations for the scalar subelement. More...
 
template<typename T >
void apply_inverse_dof_transformation (const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
 Apply inverse transformation to some data. For VectorElements, this applies the transformations for the scalar subelement. More...
 
template<typename T >
void apply_dof_transformation_to_transpose (const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
 Apply DOF transformation to some transposed data. More...
 
template<typename T >
void apply_inverse_dof_transformation_to_transpose (const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
 Apply inverse of DOF transformation to some transposed data. More...
 
template<typename T >
void apply_transpose_dof_transformation_to_transpose (const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
 Apply transpose of transformation to some transposed data. More...
 
template<typename T >
void apply_inverse_transpose_dof_transformation_to_transpose (const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
 Apply inverse transpose transformation to some transposed data. More...
 
void permute_dofs (const xtl::span< std::int32_t > &doflist, std::uint32_t cell_permutation) const
 Permute the DOFs of the element. More...
 
void unpermute_dofs (const xtl::span< std::int32_t > &doflist, std::uint32_t cell_permutation) const
 Unpermute the DOFs of the element. More...
 
std::function< void(const xtl::span< std::int32_t > &, std::uint32_t)> get_dof_permutation_function (bool inverse=false, bool scalar_element=false) const
 Return a function that applies DOF permutation to some data. More...
 

Detailed Description

Finite Element, containing the dof layout on a reference element, and various methods for evaluating and transforming the basis.

Constructor & Destructor Documentation

◆ FiniteElement() [1/2]

FiniteElement ( const ufcx_finite_element &  e)
explicit

Create finite element from UFC finite element.

Parameters
[in]eUFC finite element

◆ FiniteElement() [2/2]

FiniteElement ( const basix::FiniteElement &  element,
int  bs 
)

Create finite element from a Basix finite element.

Parameters
[in]elementBasix finite element
[in]bsThe block size

Member Function Documentation

◆ apply_dof_transformation()

void apply_dof_transformation ( const xtl::span< T > &  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply DOF transformation to some data.

Parameters
[in,out]dataThe data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_dof_transformation_to_transpose()

void apply_dof_transformation_to_transpose ( const xtl::span< T > &  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply DOF transformation to some transposed data.

Parameters
[in,out]dataThe data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_inverse_dof_transformation()

void apply_inverse_dof_transformation ( const xtl::span< T > &  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply inverse transformation to some data. For VectorElements, this applies the transformations for the scalar subelement.

Parameters
[in,out]dataThe data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_inverse_dof_transformation_to_transpose()

void apply_inverse_dof_transformation_to_transpose ( const xtl::span< T > &  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply inverse of DOF transformation to some transposed data.

Parameters
[in,out]dataThe data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_inverse_transpose_dof_transformation()

void apply_inverse_transpose_dof_transformation ( const xtl::span< T > &  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply inverse transpose transformation to some data. For VectorElements, this applies the transformations for the scalar subelement.

Parameters
[in,out]dataThe data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_inverse_transpose_dof_transformation_to_transpose()

void apply_inverse_transpose_dof_transformation_to_transpose ( const xtl::span< T > &  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply inverse transpose transformation to some transposed data.

Parameters
[in,out]dataThe data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_transpose_dof_transformation()

void apply_transpose_dof_transformation ( const xtl::span< T > &  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply transpose transformation to some data. For VectorElements, this applies the transformations for the scalar subelement.

Parameters
[in,out]dataThe data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_transpose_dof_transformation_to_transpose()

void apply_transpose_dof_transformation_to_transpose ( const xtl::span< T > &  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply transpose of transformation to some transposed data.

Parameters
[in,out]dataThe data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ block_size()

int block_size ( ) const
noexcept

Block size of the finite element function space. For VectorElements and TensorElements, this is the number of DOFs colocated at each DOF point. For other elements, this is always 1.

Returns
Block size of the finite element space

◆ cell_shape()

mesh::CellType cell_shape ( ) const
noexcept

Cell shape.

Returns
Element cell shape

◆ create_interpolation_operator()

xt::xtensor< double, 2 > create_interpolation_operator ( const FiniteElement from) const

Create a matrix that maps degrees of freedom from one element to this element (interpolation).

Parameters
[in]fromThe element to interpolate from
Returns
Matrix operator that maps the from degrees-of-freedom to the degrees-of-freedom of this element. Shape is (num_dofs of this element, num_dofs of from).
Note
The two elements must use the same mapping between the reference and physical cells
Does not support mixed elements

◆ family()

std::string family ( ) const
noexcept

The finite element family.

Returns
The string of the finite element family

◆ get_dof_permutation_function()

std::function< void(const xtl::span< std::int32_t > &, std::uint32_t)> get_dof_permutation_function ( bool  inverse = false,
bool  scalar_element = false 
) const

Return a function that applies DOF permutation to some data.

The returned function will take three inputs:

  • [in,out] doflist The numbers of the DOFs, a span of length num_dofs
  • [in] cell_permutation Permutation data for the cell
  • [in] block_size The block_size of the input data
Parameters
[in]inverseIndicates whether the inverse transformations should be returned
[in]scalar_elementIndicated whether the scalar transformations should be returned for a vector element

◆ get_dof_transformation_function()

std::function<void(const xtl::span<T>&, const xtl::span<const std::uint32_t>&, std::int32_t, int)> get_dof_transformation_function ( bool  inverse = false,
bool  transpose = false,
bool  scalar_element = false 
) const
inline

Return a function that applies DOF transformation to some data.

The returned function will take four inputs:

  • [in,out] data The data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
  • [in] cell_info Permutation data for the cell. The size of this is num_cells. For elements where no transformations are required, an empty span can be passed in.
  • [in] cell The cell number
  • [in] block_size The block_size of the input data
Parameters
[in]inverseIndicates whether the inverse transformations should be returned
[in]transposeIndicates whether the transpose transformations should be returned
[in]scalar_elementIndicates whether the scalar transformations should be returned for a vector element

◆ get_dof_transformation_to_transpose_function()

std::function<void(const xtl::span<T>&, const xtl::span<const std::uint32_t>&, std::int32_t, int)> get_dof_transformation_to_transpose_function ( bool  inverse = false,
bool  transpose = false,
bool  scalar_element = false 
) const
inline

Return a function that applies DOF transformation to some transposed data.

The returned function will take three inputs:

  • [in,out] data The data to be transformed. This data is flattened with row-major layout, shape=(num_dofs, block_size)
  • [in] cell_info Permutation data for the cell. The size of this is num_cells. For elements where no transformations are required, an empty span can be passed in.
  • [in] cell The cell number
  • [in] block_size The block_size of the input data
Parameters
[in]inverseIndicates whether the inverse transformations should be returned
[in]transposeIndicates whether the transpose transformations should be returned
[in]scalar_elementIndicated whether the scalar transformations should be returned for a vector element

◆ interpolation_ident()

bool interpolation_ident ( ) const
noexcept

Check if interpolation into the finite element space is an identity operation given the evaluation on an expression at specific points, i.e. the degree-of-freedom are equal to point evaluations. The function will return true for Lagrange elements.

Returns
True if interpolation is an identity operation

◆ interpolation_operator()

const xt::xtensor< double, 2 > & interpolation_operator ( ) const

Interpolation operator (matrix) Pi that maps a function evaluated at the points provided by FiniteElement::interpolation_points to the element degrees of freedom, i.e. dofs = Pi f_x. See the Basix documentation for basix::FiniteElement::interpolation_matrix for how the data in f_x should be ordered.

Returns
The interpolation operator Pi. Shape is (num_dofs, num_points*value_size)

◆ interpolation_points()

const xt::xtensor< double, 2 > & interpolation_points ( ) const

Points on the reference cell at which an expression need to be evaluated in order to interpolate the expression in the finite element space. For Lagrange elements the points will just be the nodal positions. For other elements the points will typically be the quadrature points used to evaluate moment degrees of freedom.

Returns
Points on the reference cell. Shape is (num_points, tdim).

◆ is_mixed()

bool is_mixed ( ) const
noexcept

Check if element is a mixed element, i.e. composed of two or more elements of different types. A block element, e.g. a Lagrange element with block size > 1 is not considered mixed.

Returns
True is element is mixed.

◆ map_fn()

std::function<void(O&, const P&, const Q&, double, const R&)> 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 pushed-forward (pulled-back) data (ndim==1)
PThe type that hold the data to be pulled back (pushed forwarded) (ndim==1)
QThe type that holds the Jacobian (inverse Jacobian) matrix (ndim==2)
RThe type that holds the inverse Jacobian (Jacobian) matrix (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 passed 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 oif the Jacobian matrix of the map, shape=(tdim, gdim)
  • detJ_inv [in] 1/det(J)
  • J [in] The Jacobian matrix, shape=(gdim, tdim)

◆ map_ident()

bool map_ident ( ) const
noexcept

Check if the push forward/pull back map from the values on reference to the values on a physical cell for this element is the identity map.

Returns
True if the map is the identity

◆ needs_dof_permutations()

bool needs_dof_permutations ( ) const
noexcept

Check if DOF permutations are needed for this element.

DOF permutations will be needed for elements which might not be continuous when two neighbouring cells disagree on the orientation of a shared subentity, and when this can be corrected for by permuting the DOF numbering in the dofmap.

For example, higher order Lagrange elements will need DOF permutations, as the arrangement of DOFs on a shared subentity may be different from the point of view of neighbouring cells, and this can be corrected for by permuting the DOF numbers on each cell.

Returns
True if DOF transformations are required

◆ needs_dof_transformations()

bool needs_dof_transformations ( ) const
noexcept

Check if DOF transformations are needed for this element.

DOF transformations will be needed for elements which might not be continuous when two neighbouring cells disagree on the orientation of a shared subentity, and when this cannot be corrected for by permuting the DOF numbering in the dofmap.

For example, Raviart-Thomas elements will need DOF transformations, as the neighbouring cells may disagree on the orientation of a basis function, and this orientation cannot be corrected for by permuting the DOF numbers on each cell.

Returns
True if DOF transformations are required

◆ num_sub_elements()

int num_sub_elements ( ) const
noexcept

Get the number of sub elements (for a mixed or blocked element)

Returns
The number of sub elements

◆ operator!=()

bool operator!= ( const FiniteElement e) const

Check if two elements are not equivalent.

Returns
True is the two elements are not the same
Note
Equality can be checked only for non-mixed elements. For a mixed element, this function will raise an exception.

◆ operator==()

bool operator== ( const FiniteElement e) const

Check if two elements are equivalent.

Returns
True is the two elements are the same
Note
Equality can be checked only for non-mixed elements. For a mixed element, this function will raise an exception.

◆ permute_dofs()

void permute_dofs ( const xtl::span< std::int32_t > &  doflist,
std::uint32_t  cell_permutation 
) const

Permute the DOFs of the element.

Parameters
[in,out]doflistThe numbers of the DOFs, a span of length num_dofs
[in]cell_permutationPermutation data for the cell

◆ reference_value_size()

int reference_value_size ( ) const

The value size, e.g. 1 for a scalar function, 2 for a 2D vector, 9 for a second-order tensor in 3D, for the reference element.

Returns
The value size for the reference element

◆ signature()

std::string signature ( ) const
noexcept

String identifying the finite element.

Returns
Element signature
Warning
The function is provided for convenience, but it should not be relied upon for determining the element type. Use other functions, commonly returning enums, to determine element properties.

◆ space_dimension()

int space_dimension ( ) const
noexcept

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

Returns
Dimension of the finite element space

◆ tabulate()

void tabulate ( xt::xtensor< double, 4 > &  values,
const xt::xtensor< double, 2 > &  X,
int  order 
) const

Evaluate all derivatives of the basis functions up to given order at given points in reference cell.

Parameters
[in,out]valuesFour dimensional xtensor that will be filled with the tabulated values. Should be of shape {num_derivatives, num_points, num_dofs, reference_value_size}
[in]XTwo dimensional xtensor of shape [num_points, geometric dimension] containing the points at the reference element
[in]orderThe number of derivatives (up to and including this order) to tabulate for

◆ unpermute_dofs()

void unpermute_dofs ( const xtl::span< std::int32_t > &  doflist,
std::uint32_t  cell_permutation 
) const

Unpermute the DOFs of the element.

Parameters
[in,out]doflistThe numbers of the DOFs, a span of length num_dofs
[in]cell_permutationPermutation data for the cell

◆ value_size()

int value_size ( ) const

The value size, e.g. 1 for a scalar function, 2 for a 2D vector, 9 for a second-order tensor in 3D.

Note
The return value of this function is equal to std::accumulate(value_shape().begin(), value_shape().end(), 1, std::multiplies<int>()).
Returns
The value size

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