Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/df/d27/classdolfinx_1_1fem_1_1FiniteElement.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
Public Member Functions | List of all members
FiniteElement< T > Class Template 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.
 
 FiniteElement (const basix::FiniteElement< T > &element, const std::vector< std::size_t > &value_shape)
 Create finite element from a Basix finite element.
 
 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.
 
bool operator!= (const FiniteElement &e) const
 Check if two elements are not equivalent.
 
std::string signature () const noexcept
 String identifying the finite element.
 
mesh::CellType cell_shape () const noexcept
 Cell shape.
 
int space_dimension () const noexcept
 Dimension of the finite element function space (the number of degrees-of-freedom for the element)
 
int block_size () const noexcept
 Block size of the finite element function space. For BlockedElements, this is the number of DOFs colocated at each DOF point. For other elements, this is always 1.
 
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.
 
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.
 
std::span< const std::size_t > 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.
 
void tabulate (std::span< T > values, std::span< const T > X, std::array< std::size_t, 2 > shape, int order) const
 Evaluate derivatives of the basis functions up to given order at points in the reference cell.
 
std::pair< std::vector< T >, std::array< std::size_t, 4 > > tabulate (std::span< const T > X, std::array< std::size_t, 2 > shape, int order) const
 Evaluate all derivatives of the basis functions up to given order at given points in reference cell.
 
int num_sub_elements () const noexcept
 Number of sub elements (for a mixed or blocked element)
 
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.
 
const std::vector< std::shared_ptr< const FiniteElement< T > > > & sub_elements () const noexcept
 Subelements (if any)
 
std::shared_ptr< const FiniteElement< T > > extract_sub_element (const std::vector< int > &component) const
 Extract sub finite element for component.
 
const basix::FiniteElement< T > & basix_element () const
 Return underlying basix element (if it exists)
 
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.
 
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.
 
std::pair< std::vector< T >, std::array< std::size_t, 2 > > interpolation_points () const
 Points on the reference cell at which an expression needs to be evaluated in order to interpolate the expression in the finite element space.
 
std::pair< std::vector< T >, std::array< std::size_t, 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.
 
std::pair< std::vector< T >, std::array< std::size_t, 2 > > create_interpolation_operator (const FiniteElement &from) const
 Create a matrix that maps degrees of freedom from one element to this element (interpolation).
 
bool needs_dof_transformations () const noexcept
 Check if DOF transformations are needed for this element.
 
bool needs_dof_permutations () const noexcept
 Check if DOF permutations are needed for this element.
 
template<typename U >
std::function< void(const std::span< U > &, const std::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.
 
template<typename U >
std::function< void(const std::span< U > &, const std::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.
 
template<typename U >
void apply_dof_transformation (const std::span< U > &data, std::uint32_t cell_permutation, int block_size) const
 Apply DOF transformation to some data.
 
template<typename U >
void apply_inverse_transpose_dof_transformation (const std::span< U > &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.
 
template<typename U >
void apply_transpose_dof_transformation (const std::span< U > &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.
 
template<typename U >
void apply_inverse_dof_transformation (const std::span< U > &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.
 
template<typename U >
void apply_dof_transformation_to_transpose (const std::span< U > &data, std::uint32_t cell_permutation, int block_size) const
 Apply DOF transformation to some transposed data.
 
template<typename U >
void apply_inverse_dof_transformation_to_transpose (const std::span< U > &data, std::uint32_t cell_permutation, int block_size) const
 Apply inverse of DOF transformation to some transposed data.
 
template<typename U >
void apply_transpose_dof_transformation_to_transpose (const std::span< U > &data, std::uint32_t cell_permutation, int block_size) const
 Apply transpose of transformation to some transposed data.
 
template<typename U >
void apply_inverse_transpose_dof_transformation_to_transpose (const std::span< U > &data, std::uint32_t cell_permutation, int block_size) const
 Apply inverse transpose transformation to some transposed data.
 
void permute_dofs (const std::span< std::int32_t > &doflist, std::uint32_t cell_permutation) const
 Permute the DOFs of the element.
 
void unpermute_dofs (const std::span< std::int32_t > &doflist, std::uint32_t cell_permutation) const
 Unpermute the DOFs of the element.
 
std::function< void(const std::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.
 

Detailed Description

template<std::floating_point T>
class dolfinx::fem::FiniteElement< T >

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]

template<std::floating_point T>
FiniteElement ( const ufcx_finite_element &  e)
explicit

Create finite element from UFC finite element.

Parameters
[in]eUFC finite element.

◆ FiniteElement() [2/2]

template<std::floating_point T>
FiniteElement ( const basix::FiniteElement< T > &  element,
const std::vector< std::size_t > &  value_shape 
)

Create finite element from a Basix finite element.

Parameters
[in]elementBasix finite element
[in]value_shapeValue shape for 'blocked' elements, e.g. vector-valued Lagrange elements where each component for the vector field is a Lagrange element. For example, a vector-valued element in 3D will have value_shape equal to {3}, and for a second-order tensor element in 2D value_shape equal to {2, 2}.

Member Function Documentation

◆ apply_dof_transformation()

template<std::floating_point T>
template<typename U >
void apply_dof_transformation ( const std::span< U > &  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()

template<std::floating_point T>
template<typename U >
void apply_dof_transformation_to_transpose ( const std::span< U > &  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()

template<std::floating_point T>
template<typename U >
void apply_inverse_dof_transformation ( const std::span< U > &  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()

template<std::floating_point T>
template<typename U >
void apply_inverse_dof_transformation_to_transpose ( const std::span< U > &  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()

template<std::floating_point T>
template<typename U >
void apply_inverse_transpose_dof_transformation ( const std::span< U > &  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()

template<std::floating_point T>
template<typename U >
void apply_inverse_transpose_dof_transformation_to_transpose ( const std::span< U > &  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()

template<std::floating_point T>
template<typename U >
void apply_transpose_dof_transformation ( const std::span< U > &  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()

template<std::floating_point T>
template<typename U >
void apply_transpose_dof_transformation_to_transpose ( const std::span< U > &  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()

template<std::floating_point T>
int block_size ( ) const
noexcept

Block size of the finite element function space. For BlockedElements, 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()

template<std::floating_point T>
mesh::CellType cell_shape ( ) const
noexcept

Cell shape.

Returns
Element cell shape

◆ create_interpolation_operator()

template<std::floating_point T>
std::pair< std::vector< T >, std::array< std::size_t, 2 > > create_interpolation_operator ( const FiniteElement< T > &  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. The (0) matrix data (row-major storage) and (1) the shape (num_dofs of this element, num_dofs of from) are returned.
Precondition
The two elements must use the same mapping between the reference and physical cells
Note
Does not support mixed elements

◆ family()

template<std::floating_point T>
std::string family ( ) const
noexcept

The finite element family.

Returns
The string of the finite element family

◆ get_dof_permutation_function()

template<std::floating_point T>
std::function< void(const std::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()

template<std::floating_point T>
template<typename U >
std::function< void(const std::span< U > &, const std::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()

template<std::floating_point T>
template<typename U >
std::function< void(const std::span< U > &, const std::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()

template<std::floating_point T>
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()

template<std::floating_point T>
std::pair< std::vector< T >, std::array< std::size_t, 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, returning the data for Pi (row-major storage) and the shape (num_dofs, num_points * value_size)

◆ interpolation_points()

template<std::floating_point T>
std::pair< std::vector< T >, std::array< std::size_t, 2 > > interpolation_points ( ) const

Points on the reference cell at which an expression needs 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
Interpolation point coordinates on the reference cell, returning the (0) coordinates data (row-major) storage and (1) the shape (num_points, tdim).

◆ is_mixed()

template<std::floating_point T>
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_ident()

template<std::floating_point T>
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()

template<std::floating_point T>
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 sub-entity 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()

template<std::floating_point T>
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 sub-entity, 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()

template<std::floating_point T>
int num_sub_elements ( ) const
noexcept

Number of sub elements (for a mixed or blocked element)

Returns
The number of sub elements

◆ operator!=()

template<std::floating_point T>
bool operator!= ( const FiniteElement< T > &  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==()

template<std::floating_point T>
bool operator== ( const FiniteElement< T > &  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()

template<std::floating_point T>
void permute_dofs ( const std::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()

template<std::floating_point T>
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()

template<std::floating_point T>
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()

template<std::floating_point T>
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() [1/2]

template<std::floating_point T>
std::pair< std::vector< T >, std::array< std::size_t, 4 > > tabulate ( std::span< const T >  X,
std::array< std::size_t, 2 >  shape,
int  order 
) const

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

Parameters
[in]XThe reference coordinates at which to evaluate the basis functions. Shape is (num_points, topological dimension) (row-major storage)
[in]shapeThe shape of X
[in]orderThe number of derivatives (up to and including this order) to tabulate for
Returns
Basis function values and array shape (row-major storage)

◆ tabulate() [2/2]

template<std::floating_point T>
void tabulate ( std::span< T >  values,
std::span< const T >  X,
std::array< std::size_t, 2 >  shape,
int  order 
) const

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

Parameters
[in,out]valuesArray that will be filled with the tabulated basis values. Must have shape (num_derivatives, num_points, num_dofs, reference_value_size) (row-major storage)
[in]XThe reference coordinates at which to evaluate the basis functions. Shape is (num_points, topological dimension) (row-major storage)
[in]shapeThe shape of X
[in]orderThe number of derivatives (up to and including this order) to tabulate for

◆ unpermute_dofs()

template<std::floating_point T>
void unpermute_dofs ( const std::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()

template<std::floating_point T>
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 equivalent to std::accumulate(value_shape().begin(), value_shape().end(), 1, std::multiplies{}).
Returns
The value size

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