Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.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 evaluate_reference_basis (xt::xtensor< double, 3 > &values, const xt::xtensor< double, 2 > &X) const
 Evaluate all basis functions at given points in reference cell.
 
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
 Evaluate all basis function derivatives of given order at given points in reference cell. More...
 
void transform_reference_basis_derivatives (std::vector< double > &values, std::size_t order, const std::vector< double > &reference_values, const array2d< double > &X, const std::vector< double > &J, const xtl::span< const double > &detJ, const std::vector< double > &K) const
 Push basis function (derivatives) 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_permutation_data () const noexcept
 
template<typename T >
void apply_dof_transformation (xtl::span< T > data, std::uint32_t cell_permutation, int block_size) const
 Apply permutation to some data. More...
 
template<typename T >
void apply_inverse_transpose_dof_transformation (xtl::span< T > data, std::uint32_t cell_permutation, int block_size) const
 Apply inverse transpose permutation to some data. More...
 
template<typename T >
void map_pull_back (const xt::xtensor< T, 3 > &u, const xt::xtensor< double, 3 > &J, const xtl::span< const double > &detJ, const xt::xtensor< double, 3 > &K, xt::xtensor< T, 3 > &U) const
 Pull physical data back to the reference element. This passes the inputs directly into Basix's map_pull_back function.
 

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 ( xtl::span< T >  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply permutation 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_inverse_transpose_dof_transformation()

template<typename T >
void dolfinx::fem::FiniteElement::apply_inverse_transpose_dof_transformation ( xtl::span< T >  data,
std::uint32_t  cell_permutation,
int  block_size 
) const
inline

Apply inverse transpose permutation 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

◆ 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

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

◆ needs_permutation_data()

bool FiniteElement::needs_permutation_data ( ) const
noexcept
Todo:
Expand on when permutation data might be required

Check if cell permutation data is required for this element

Returns
True if cell permutation data is 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

◆ 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

◆ transform_reference_basis()

void FiniteElement::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

Evaluate all basis function derivatives of given order at given points in reference cell.

Push basis functions forward to physical element

◆ 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: