11#include <basix/finite-element.h>
14#include <dolfinx/mesh/cell_types.h>
21struct ufcx_finite_element;
38template <std::
floating_po
int T>
130 std::array<std::
size_t, 2> shape,
int order) const;
141 std::pair<std::vector<
geometry_type>, std::array<std::
size_t, 4>>
197 std::pair<std::vector<
geometry_type>, std::array<std::
size_t, 2>>
209 std::pair<std::vector<
geometry_type>, std::array<std::
size_t, 2>>
224 std::pair<std::vector<
geometry_type>, std::array<std::
size_t, 2>>
297 template <typename U>
298 std::function<
void(std::span<U>, std::span<const std::uint32_t>, std::int32_t,
305 return [](std::span<U>, std::span<const std::uint32_t>, std::int32_t, int)
311 if (!_sub_elements.empty())
316 std::vector<std::function<void(
317 std::span<U>, std::span<const std::uint32_t>, std::int32_t,
int)>>
318 sub_element_functions;
319 std::vector<int> dims;
320 for (std::size_t i = 0; i < _sub_elements.size(); ++i)
322 sub_element_functions.push_back(
323 _sub_elements[i]->
template dof_transformation_fn<U>(ttype));
327 return [dims, sub_element_functions](
328 std::span<U> data, std::span<const std::uint32_t> cell_info,
331 std::size_t offset = 0;
332 for (std::size_t e = 0; e < sub_element_functions.size(); ++e)
334 const std::size_t width = dims[e] *
block_size;
335 sub_element_functions[e](data.subspan(offset, width), cell_info,
341 else if (!scalar_element)
344 std::function<void(std::span<U>, std::span<const std::uint32_t>,
347 = _sub_elements[0]->template dof_transformation_fn<U>(ttype);
349 return [ebs, sub_function](std::span<U> data,
350 std::span<const std::uint32_t> cell_info,
351 std::int32_t
cell,
int data_block_size)
352 { sub_function(data, cell_info,
cell, ebs * data_block_size); };
359 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
363 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
367 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
371 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
375 throw std::runtime_error(
"Unknown transformation type");
402 template <
typename U>
403 std::function<void(std::span<U>, std::span<const std::uint32_t>, std::int32_t,
406 bool scalar_element =
false)
const
411 return [](std::span<U>, std::span<const std::uint32_t>, std::int32_t, int)
416 else if (_sub_elements.size() != 0)
421 std::vector<std::function<void(
422 std::span<U>, std::span<const std::uint32_t>, std::int32_t,
int)>>
423 sub_element_functions;
424 for (std::size_t i = 0; i < _sub_elements.size(); ++i)
426 sub_element_functions.push_back(
427 _sub_elements[i]->
template dof_transformation_right_fn<U>(ttype));
430 return [
this, sub_element_functions](
431 std::span<U> data, std::span<const std::uint32_t> cell_info,
434 std::size_t offset = 0;
435 for (std::size_t e = 0; e < sub_element_functions.size(); ++e)
437 sub_element_functions[e](data.subspan(offset, data.size() - offset),
439 offset += _sub_elements[e]->space_dimension();
443 else if (!scalar_element)
450 std::function<void(std::span<U>, std::span<const std::uint32_t>,
453 = _sub_elements[0]->template dof_transformation_fn<U>(ttype);
454 return [
this, sub_function](std::span<U> data,
455 std::span<const std::uint32_t> cell_info,
456 std::int32_t
cell,
int data_block_size)
459 const std::size_t dof_count = data.size() / data_block_size;
460 for (
int block = 0; block < data_block_size; ++block)
462 sub_function(data.subspan(block * dof_count, dof_count), cell_info,
472 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
473 std::int32_t
cell,
int n)
476 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
477 std::int32_t
cell,
int n)
480 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
481 std::int32_t
cell,
int n)
484 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
485 std::int32_t
cell,
int n)
488 throw std::runtime_error(
"Unknown transformation type");
528 template <
typename U>
529 void T_apply(std::span<U> data, std::uint32_t cell_permutation,
int n)
const
532 _element->T_apply(data, n, cell_permutation);
544 template <
typename U>
549 _element->Tt_inv_apply(data, n, cell_permutation);
561 template <
typename U>
562 void Tt_apply(std::span<U> data, std::uint32_t cell_permutation,
int n)
const
565 _element->Tt_apply(data, n, cell_permutation);
576 template <
typename U>
577 void Tinv_apply(std::span<U> data, std::uint32_t cell_permutation,
581 _element->Tinv_apply(data, n, cell_permutation);
592 template <
typename U>
597 _element->T_apply_right(data, n, cell_permutation);
609 template <
typename U>
614 _element->Tinv_apply_right(data, n, cell_permutation);
626 template <
typename U>
631 _element->Tt_apply_right(data, n, cell_permutation);
643 template <
typename U>
648 _element->Tt_inv_apply_right(data, n, cell_permutation);
667 void permute(std::span<std::int32_t> doflist,
668 std::uint32_t cell_permutation)
const;
687 std::uint32_t cell_permutation)
const;
703 std::function<void(std::span<std::int32_t>, std::uint32_t)>
707 std::string _signature;
714 std::vector<std::shared_ptr<const FiniteElement<geometry_type>>>
718 std::vector<std::size_t> _reference_value_shape;
731 bool _needs_dof_permutations;
732 bool _needs_dof_transformations;
735 std::unique_ptr<basix::FiniteElement<geometry_type>> _element;
739 std::pair<std::vector<geometry_type>, std::array<std::size_t, 2>> _points;
Definition CoordinateElement.h:26
Model of a finite element.
Definition FiniteElement.h:40
bool symmetric() const
Does the element represent a symmetric 2-tensor?
Definition FiniteElement.cpp:382
FiniteElement & operator=(const FiniteElement &element)=delete
Copy assignment.
const std::vector< std::shared_ptr< const FiniteElement< geometry_type > > > & sub_elements() const noexcept
Get subelements (if any)
Definition FiniteElement.cpp:432
void Tinv_apply(std::span< U > data, std::uint32_t cell_permutation, int n) const
Apply the inverse of the operator applied by T_apply().
Definition FiniteElement.h:577
std::pair< std::vector< geometry_type >, 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).
Definition FiniteElement.cpp:522
void T_apply(std::span< U > data, std::uint32_t cell_permutation, int n) const
Transform basis functions from the reference element ordering and orientation to the globally consist...
Definition FiniteElement.h:529
basix::maps::type map_type() const
Get the map type used by the element.
Definition FiniteElement.cpp:461
int space_dimension() const noexcept
Definition FiniteElement.cpp:370
void tabulate(std::span< geometry_type > values, std::span< const geometry_type > 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.
Definition FiniteElement.cpp:401
void Tt_inv_apply_right(std::span< U > data, std::uint32_t cell_permutation, int n) const
Right(post)-apply the transpose inverse of the operator applied by T_apply().
Definition FiniteElement.h:644
FiniteElement & operator=(FiniteElement &&element)=default
Move assignment.
bool needs_dof_transformations() const noexcept
Check if DOF transformations are needed for this element.
Definition FiniteElement.cpp:565
FiniteElement(const FiniteElement &element)=delete
Copy constructor.
FiniteElement(const ufcx_finite_element &e)
Create finite element from UFC finite element.
Definition FiniteElement.cpp:124
void Tt_inv_apply(std::span< U > data, std::uint32_t cell_permutation, int n) const
Apply the inverse transpose of the operator applied by T_apply().
Definition FiniteElement.h:545
std::pair< std::vector< geometry_type >, 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...
Definition FiniteElement.cpp:493
void permute_inv(std::span< std::int32_t > doflist, std::uint32_t cell_permutation) const
Perform the inverse of the operation applied by permute().
Definition FiniteElement.cpp:584
std::shared_ptr< const FiniteElement< geometry_type > > extract_sub_element(const std::vector< int > &component) const
Extract sub finite element for component.
Definition FiniteElement.cpp:439
void T_apply_right(std::span< U > data, std::uint32_t cell_permutation, int n) const
Right(post)-apply the operator applied by T_apply().
Definition FiniteElement.h:593
std::string signature() const noexcept
Definition FiniteElement.cpp:358
bool needs_dof_permutations() const noexcept
Check if DOF permutations are needed for this element.
Definition FiniteElement.cpp:571
bool is_mixed() const noexcept
Check if element is a mixed element.
Definition FiniteElement.cpp:425
int reference_value_size() const
Definition FiniteElement.cpp:388
std::function< void(std::span< U >, std::span< const std::uint32_t >, std::int32_t, int)> dof_transformation_fn(doftransform ttype, bool scalar_element=false) const
Return a function that applies a DOF transformation operator to some data (see T_apply()).
Definition FiniteElement.h:300
bool operator!=(const FiniteElement &e) const
Definition FiniteElement.cpp:352
~FiniteElement()=default
Destructor.
T geometry_type
Geometry type of the Mesh that the FunctionSpace is defined on.
Definition FiniteElement.h:43
const basix::FiniteElement< geometry_type > & basix_element() const
Return underlying Basix element (if it exists).
Definition FiniteElement.cpp:449
std::function< void(std::span< std::int32_t >, std::uint32_t)> dof_permutation_fn(bool inverse=false, bool scalar_element=false) const
Return a function that applies a degree-of-freedom permutation to some data.
FiniteElement(FiniteElement &&element)=default
Move constructor.
int num_sub_elements() const noexcept
Number of sub elements (for a mixed or blocked element).
Definition FiniteElement.cpp:419
void Tt_apply(std::span< U > data, std::uint32_t cell_permutation, int n) const
Apply the transpose of the operator applied by T_apply().
Definition FiniteElement.h:562
std::pair< std::vector< geometry_type >, std::array< std::size_t, 2 > > interpolation_operator() const
Definition FiniteElement.cpp:509
int block_size() const noexcept
Definition FiniteElement.cpp:395
std::span< const std::size_t > reference_value_shape() const
The reference value shape.
Definition FiniteElement.cpp:376
void Tt_apply_right(std::span< U > data, std::uint32_t cell_permutation, int n) const
Right(post)-apply the transpose of the operator applied by T_apply().
Definition FiniteElement.h:627
std::function< void(std::span< U >, std::span< const std::uint32_t >, std::int32_t, int)> dof_transformation_right_fn(doftransform ttype, bool scalar_element=false) const
Return a function that applies DOF transformation to some transposed data (see T_apply_right()).
Definition FiniteElement.h:405
void permute(std::span< std::int32_t > doflist, std::uint32_t cell_permutation) const
Permute indices associated with degree-of-freedoms on the reference element ordering to the globally ...
Definition FiniteElement.cpp:577
bool interpolation_ident() const noexcept
Definition FiniteElement.cpp:483
void Tinv_apply_right(std::span< U > data, std::uint32_t cell_permutation, int n) const
Right(post)-apply the inverse of the operator applied by T_apply().
Definition FiniteElement.h:610
bool operator==(const FiniteElement &e) const
Definition FiniteElement.cpp:340
mesh::CellType cell_shape() const noexcept
Definition FiniteElement.cpp:364
bool map_ident() const noexcept
Definition FiniteElement.cpp:473
Finite element method functionality.
Definition assemble_matrix_impl.h:26
doftransform
DOF transformation type.
Definition FiniteElement.h:27
@ inverse_transpose
Transpose inverse.
CellType
Cell type identifier.
Definition cell_types.h:22