Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.3.0/v0.9.0/cpp
DOLFINx  0.3.0
DOLFINx C++ interface
Public Member Functions | List of all members
dolfinx::fem::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 ufc_finite_element &ufc_element)
 Create finite element from UFC 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.
 
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. 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 noexcept
 The value size, e.g. 1 for a scalar function, 2 for a 2D vector. More...
 
int reference_value_size () const noexcept
 The value size, e.g. 1 for a scalar function, 2 for a 2D vector for the reference element. More...
 
int value_rank () const noexcept
 Rank of the value space. More...
 
int value_dimension (int i) const
 Return the dimension of the value space for axis i.
 
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...
 
void transform_reference_basis (xt::xtensor< double, 3 > &values, const xt::xtensor< double, 3 > &reference_values, const xt::xtensor< double, 3 > &J, const xtl::span< const double > &detJ, const xt::xtensor< double, 3 > &K) const
 Push basis functions forward to physical element.
 
int num_sub_elements () const noexcept
 Get the number of sub elements (for a mixed element) More...
 
const std::vector< std::shared_ptr< const FiniteElement > > & sub_elements () const noexcept
 Subelements (if any)
 
std::size_t hash () const noexcept
 Return simple hash of the signature string.
 
std::shared_ptr< const FiniteElementextract_sub_element (const std::vector< int > &component) const
 Extract sub finite element for component.
 
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...
 
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...
 
template<typename T >
constexpr void interpolate (const xt::xtensor< T, 2 > &values, xtl::span< T > dofs) const
 
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 tranposed 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...
 
template<typename O , typename P , typename Q , typename T , typename S >
void map_pull_back (const O &u, const P &J, const Q &detJ, const T &K, S &&U) const
 Pull back data from the physical element to the reference element. It can process batches of points that share the same geometric map. 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 transformation 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()

FiniteElement::FiniteElement ( const ufc_finite_element &  ufc_element)
explicit

Create finite element from UFC finite element.

Parameters
[in]ufc_elementUFC finite element

Member Function Documentation

◆ apply_dof_transformation()

template<typename T >
void dolfinx::fem::FiniteElement::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
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_dof_transformation_to_transpose()

template<typename T >
void dolfinx::fem::FiniteElement::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 tranposed data.

Parameters
[in,out]dataThe data to be transformed
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_inverse_dof_transformation()

template<typename T >
void dolfinx::fem::FiniteElement::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
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_inverse_dof_transformation_to_transpose()

template<typename T >
void dolfinx::fem::FiniteElement::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
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_inverse_transpose_dof_transformation()

template<typename T >
void dolfinx::fem::FiniteElement::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
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_inverse_transpose_dof_transformation_to_transpose()

template<typename T >
void dolfinx::fem::FiniteElement::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
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_transpose_dof_transformation()

template<typename T >
void dolfinx::fem::FiniteElement::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
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ apply_transpose_dof_transformation_to_transpose()

template<typename T >
void dolfinx::fem::FiniteElement::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
[in]cell_permutationPermutation data for the cell
[in]block_sizeThe block_size of the input data

◆ block_size()

int FiniteElement::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 FiniteElement::cell_shape ( ) const
noexcept

Cell shape.

Returns
Element cell shape

◆ family()

std::string FiniteElement::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)> FiniteElement::get_dof_permutation_function ( bool  inverse = false,
bool  scalar_element = false 
) const

Return a function that applies DOF transformation to some data.

The returned function will take three inputs:

  • [in,out] data The data to be transformed
  • [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<typename T >
std::function<void(const xtl::span<T>&, const xtl::span<const std::uint32_t>&, std::int32_t, int)> dolfinx::fem::FiniteElement::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 three inputs:

  • [in,out] data The data to be transformed
  • [in] cell_info Permutation data for the cell
  • [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

◆ get_dof_transformation_to_transpose_function()

template<typename T >
std::function<void(const xtl::span<T>&, const xtl::span<const std::uint32_t>&, std::int32_t, int)> dolfinx::fem::FiniteElement::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
  • [in] cell_info Permutation data for the cell
  • [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

◆ interpolate()

template<typename T >
constexpr void dolfinx::fem::FiniteElement::interpolate ( const xt::xtensor< T, 2 > &  values,
xtl::span< T >  dofs 
) const
inlineconstexpr
Todo:

Document shape/layout of values

Make the interpolating dofs in/out argument for efficiency as this function is often called from within tight loops

Consider handling block size > 1

Re-work for fields that require a pull-back, e.g. Piols mapped elements

Interpolate a function in the finite element space on a cell. Given the evaluation of the function to be interpolated at points provided by FiniteElement::interpolation_points, it evaluates the degrees of freedom for the interpolant.

Parameters
[in]valuesThe values of the function. It has shape (value_size, num_points), where num_points is the number of points given by FiniteElement::interpolation_points.
[out]dofsThe element degrees of freedom (interpolants) of the expression. The call must allocate the space. Is has

◆ interpolation_ident()

bool FiniteElement::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_points()

const xt::xtensor< double, 2 > & FiniteElement::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).

◆ map_pull_back()

template<typename O , typename P , typename Q , typename T , typename S >
void dolfinx::fem::FiniteElement::map_pull_back ( const O &  u,
const P &  J,
const Q &  detJ,
const T &  K,
S &&  U 
) const
inline

Pull back data from the physical element to the reference element. It can process batches of points that share the same geometric map.

Note
This passes the inputs directly into the Basix map_pull_back function
Parameters
[in]uData defined on the physical element. It must have dimension 3. The first index is for the geometric/map data, the second is the point index for points that share map data, and the third index is (vector) component, e.g. u[i,:,:] are points that are mapped by J[i,:,:].
[in]JThe Jacobians. It must have dimension 3. The first index is for the ith Jacobian, i.e. J[i,:,:] is the ith Jacobian.
[in]detJThe determinant of J. detJ[i] is det(J[i,:,:]). It must have dimension 1.
[in]KThe inverse of J, K[i,:,:] = J[i,:,:]^-1. It must have dimension 3.
[out]UThe input u mapped to the reference element. It must have dimension 3.

◆ needs_dof_permutations()

bool FiniteElement::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 FiniteElement::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 FiniteElement::num_sub_elements ( ) const
noexcept

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

Returns
the Number of sub elements

◆ permute_dofs()

void FiniteElement::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
[in]cell_permutationPermutation data for the cell

◆ reference_value_size()

int FiniteElement::reference_value_size ( ) const
noexcept

The value size, e.g. 1 for a scalar function, 2 for a 2D vector for the reference element.

Returns
The value size for the reference element

◆ signature()

std::string FiniteElement::signature ( ) const
noexcept

String identifying the finite element.

Returns
Element signature

◆ space_dimension()

int FiniteElement::space_dimension ( ) const
noexcept

Dimension of the finite element function space.

Returns
Dimension of the finite element space

◆ tabulate()

void FiniteElement::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 FiniteElement::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
[in]cell_permutationPermutation data for the cell

◆ value_rank()

int FiniteElement::value_rank ( ) const
noexcept

Rank of the value space.

Returns
The value rank

◆ value_size()

int FiniteElement::value_size ( ) const
noexcept

The value size, e.g. 1 for a scalar function, 2 for a 2D vector.

Returns
The value size

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