DOLFINx
0.3.0
DOLFINx C++ interface
|
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. | |
FiniteElement & | operator= (const FiniteElement &element)=delete |
Copy assignment. | |
FiniteElement & | operator= (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 FiniteElement > | extract_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... | |
Finite Element, containing the dof layout on a reference element, and various methods for evaluating and transforming the basis.
|
explicit |
Create finite element from UFC finite element.
[in] | ufc_element | UFC finite element |
|
inline |
Apply DOF transformation to some data.
[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 |
|
inline |
Apply DOF transformation to some tranposed data.
[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 |
|
inline |
Apply inverse transformation to some data. For VectorElements, this applies the transformations for the scalar subelement.
[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 |
|
inline |
Apply inverse of DOF transformation to some transposed data.
[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 |
|
inline |
Apply inverse transpose transformation to some data. For VectorElements, this applies the transformations for the scalar subelement.
[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 |
|
inline |
Apply inverse transpose transformation to some transposed data.
[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 |
|
inline |
Apply transpose transformation to some data. For VectorElements, this applies the transformations for the scalar subelement.
[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 |
|
inline |
Apply transpose of transformation to some transposed data.
[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 |
|
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.
|
noexcept |
Cell shape.
|
noexcept |
The finite element family.
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] | inverse | Indicates whether the inverse transformations should be returned |
[in] | scalar_element | Indicated whether the scalar transformations should be returned for a vector element |
|
inline |
Return a function that applies DOF transformation to some data.
The returned function will take three inputs:
[in] | inverse | Indicates whether the inverse transformations should be returned |
[in] | transpose | Indicates whether the transpose transformations should be returned |
[in] | scalar_element | Indicated whether the scalar transformations should be returned for a vector element |
|
inline |
Return a function that applies DOF transformation to some transposed data.
The returned function will take three inputs:
[in] | inverse | Indicates whether the inverse transformations should be returned |
[in] | transpose | Indicates whether the transpose transformations should be returned |
[in] | scalar_element | Indicated whether the scalar transformations should be returned for a vector element |
|
inlineconstexpr |
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.
[in] | values | The 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] | dofs | The element degrees of freedom (interpolants) of the expression. The call must allocate the space. Is has |
|
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.
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.
|
inline |
Pull back data from the physical element to the reference element. It can process batches of points that share the same geometric map.
map_pull_back
function[in] | u | Data 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] | J | The Jacobians. It must have dimension 3. The first index is for the ith Jacobian, i.e. J[i,:,:] is the ith Jacobian. |
[in] | detJ | The determinant of J. detJ[i] is det(J[i,:,:]) . It must have dimension 1. |
[in] | K | The inverse of J, K[i,:,:] = J[i,:,:]^-1 . It must have dimension 3. |
[out] | U | The input u mapped to the reference element. It must have dimension 3. |
|
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.
|
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.
|
noexcept |
Get the number of sub elements (for a mixed element)
void FiniteElement::permute_dofs | ( | const xtl::span< std::int32_t > & | doflist, |
std::uint32_t | cell_permutation | ||
) | const |
Permute the DOFs of the element.
[in,out] | doflist | The numbers of the DOFs |
[in] | cell_permutation | Permutation data for the cell |
|
noexcept |
The value size, e.g. 1 for a scalar function, 2 for a 2D vector for the reference element.
|
noexcept |
String identifying the finite element.
|
noexcept |
Dimension of the finite element function space.
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.
[in,out] | values | Four dimensional xtensor that will be filled with the tabulated values. Should be of shape {num_derivatives, num_points, num_dofs, reference_value_size} |
[in] | X | Two dimensional xtensor of shape [num_points, geometric dimension] containing the points at the reference element |
[in] | order | The number of derivatives (up to and including this order) to tabulate for. |
void FiniteElement::unpermute_dofs | ( | const xtl::span< std::int32_t > & | doflist, |
std::uint32_t | cell_permutation | ||
) | const |
Unpermute the DOFs of the element.
[in,out] | doflist | The numbers of the DOFs |
[in] | cell_permutation | Permutation data for the cell |
|
noexcept |
Rank of the value space.
|
noexcept |
The value size, e.g. 1 for a scalar function, 2 for a 2D vector.