11#include <basix/finite-element.h>
14#include <dolfinx/mesh/cell_types.h>
36template <std::
floating_po
int T>
63 std::array<std::size_t, 2> pshape, std::size_t
block_size,
121 const std::vector<std::vector<std::vector<
int>>>&
125 const std::vector<std::vector<std::vector<
int>>>&
144 std::array<std::
size_t, 2> shape,
int order) const;
155 std::pair<std::vector<
geometry_type>, std::array<std::
size_t, 4>>
211 std::pair<std::vector<
geometry_type>, std::array<std::
size_t, 2>>
223 std::pair<std::vector<
geometry_type>, std::array<std::
size_t, 2>>
238 std::pair<std::vector<
geometry_type>, std::array<std::
size_t, 2>>
311 template <typename U>
312 std::function<
void(std::span<U>, std::span<const std::uint32_t>, std::int32_t,
319 return [](std::span<U>, std::span<const std::uint32_t>, std::int32_t, int)
325 if (!_sub_elements.empty())
330 std::vector<std::function<void(
331 std::span<U>, std::span<const std::uint32_t>, std::int32_t, int)>>
333 std::vector<int> dims;
334 for (std::size_t i = 0; i < _sub_elements.size(); ++i)
336 sub_element_fns.push_back(
337 _sub_elements[i]->template dof_transformation_fn<U>(ttype));
338 dims.push_back(_sub_elements[i]->space_dimension());
341 return [dims, sub_element_fns](std::span<U> data,
342 std::span<const std::uint32_t> cell_info,
343 std::int32_t cell, int block_size)
345 std::size_t offset = 0;
346 for (std::size_t e = 0; e < sub_element_fns.size(); ++e)
348 const std::size_t width = dims[e] * block_size;
349 sub_element_fns[e](data.subspan(offset, width), cell_info, cell,
355 else if (!scalar_element)
358 std::function<void(std::span<U>, std::span<const std::uint32_t>,
362 return [ebs, sub_fn](std::span<U> data,
363 std::span<const std::uint32_t> cell_info,
364 std::int32_t
cell,
int data_block_size)
365 { sub_fn(data, cell_info,
cell, ebs * data_block_size); };
372 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
376 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
380 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
384 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
388 throw std::runtime_error(
"Unknown transformation type");
415 template <
typename U>
416 std::function<void(std::span<U>, std::span<const std::uint32_t>, std::int32_t,
419 bool scalar_element =
false)
const
421 if (!needs_dof_transformations())
424 return [](std::span<U>, std::span<const std::uint32_t>, std::int32_t, int)
429 else if (_sub_elements.size() != 0)
434 std::vector<std::function<void(
435 std::span<U>, std::span<const std::uint32_t>, std::int32_t,
int)>>
437 for (std::size_t i = 0; i < _sub_elements.size(); ++i)
439 sub_element_fns.push_back(
440 _sub_elements[i]->
template dof_transformation_right_fn<U>(ttype));
443 return [
this, sub_element_fns](std::span<U> data,
444 std::span<const std::uint32_t> cell_info,
445 std::int32_t
cell,
int block_size)
447 std::size_t offset = 0;
448 for (std::size_t e = 0; e < sub_element_fns.size(); ++e)
450 sub_element_fns[e](data.subspan(offset, data.size() - offset),
451 cell_info,
cell, block_size);
452 offset += _sub_elements[e]->space_dimension();
456 else if (!scalar_element)
463 std::function<void(std::span<U>, std::span<const std::uint32_t>,
465 sub_fn = _sub_elements[0]->template dof_transformation_fn<U>(ttype);
466 return [
this, sub_fn](std::span<U> data,
467 std::span<const std::uint32_t> cell_info,
468 std::int32_t
cell,
int data_block_size)
470 const int ebs = block_size();
471 const std::size_t dof_count = data.size() / data_block_size;
472 for (
int block = 0; block < data_block_size; ++block)
474 sub_fn(data.subspan(block * dof_count, dof_count), cell_info,
cell,
484 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
485 std::int32_t
cell,
int n)
486 { Tt_inv_apply_right(data, cell_info[
cell], n); };
488 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
489 std::int32_t
cell,
int n)
490 { Tt_apply_right(data, cell_info[
cell], n); };
492 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
493 std::int32_t
cell,
int n)
494 { Tinv_apply_right(data, cell_info[
cell], n); };
496 return [
this](std::span<U> data, std::span<const std::uint32_t> cell_info,
497 std::int32_t
cell,
int n)
498 { T_apply_right(data, cell_info[
cell], n); };
500 throw std::runtime_error(
"Unknown transformation type");
540 template <
typename U>
541 void T_apply(std::span<U> data, std::uint32_t cell_permutation,
int n)
const
544 _element->T_apply(data, n, cell_permutation);
556 template <
typename U>
561 _element->Tt_inv_apply(data, n, cell_permutation);
573 template <
typename U>
574 void Tt_apply(std::span<U> data, std::uint32_t cell_permutation,
int n)
const
577 _element->Tt_apply(data, n, cell_permutation);
588 template <
typename U>
589 void Tinv_apply(std::span<U> data, std::uint32_t cell_permutation,
593 _element->Tinv_apply(data, n, cell_permutation);
604 template <
typename U>
609 _element->T_apply_right(data, n, cell_permutation);
621 template <
typename U>
626 _element->Tinv_apply_right(data, n, cell_permutation);
638 template <
typename U>
643 _element->Tt_apply_right(data, n, cell_permutation);
655 template <
typename U>
660 _element->Tt_inv_apply_right(data, n, cell_permutation);
679 void permute(std::span<std::int32_t> doflist,
680 std::uint32_t cell_permutation)
const;
698 void permute_inv(std::span<std::int32_t> doflist,
699 std::uint32_t cell_permutation)
const;
715 std::function<void(std::span<std::int32_t>, std::uint32_t)>
719 std::string _signature;
724 std::vector<std::shared_ptr<const FiniteElement<geometry_type>>>
728 std::vector<std::size_t> _reference_value_shape;
738 std::unique_ptr<basix::FiniteElement<geometry_type>> _element;
744 bool _needs_dof_permutations;
745 bool _needs_dof_transformations;
747 std::vector<std::vector<std::vector<int>>> _entity_dofs;
748 std::vector<std::vector<std::vector<int>>> _entity_closure_dofs;
752 std::pair<std::vector<geometry_type>, std::array<std::size_t, 2>> _points;
Model of a finite element.
Definition FiniteElement.h:38
bool symmetric() const
Does the element represent a symmetric 2-tensor?
Definition FiniteElement.cpp:250
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:300
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:589
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:390
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:541
basix::maps::type map_type() const
Get the map type used by the element.
Definition FiniteElement.cpp:329
int space_dimension() const noexcept
Definition FiniteElement.cpp:223
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:269
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:656
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:433
FiniteElement(const FiniteElement &element)=delete
Copy constructor.
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:557
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:361
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:307
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:605
std::string signature() const noexcept
Definition FiniteElement.cpp:217
bool needs_dof_permutations() const noexcept
Check if DOF permutations are needed for this element.
Definition FiniteElement.cpp:439
const std::vector< std::vector< std::vector< int > > > & entity_dofs() const noexcept
The local DOFs associated with each subentity of the cell.
Definition FiniteElement.cpp:237
bool is_mixed() const noexcept
Check if element is a mixed element.
Definition FiniteElement.cpp:293
int reference_value_size() const
Definition FiniteElement.cpp:256
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:314
bool operator!=(const FiniteElement &e) const
Definition FiniteElement.cpp:211
~FiniteElement()=default
Destructor.
T geometry_type
Geometry type of the Mesh that the FunctionSpace is defined on.
Definition FiniteElement.h:41
const basix::FiniteElement< geometry_type > & basix_element() const
Return underlying Basix element (if it exists).
Definition FiniteElement.cpp:317
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(const basix::FiniteElement< geometry_type > &element, std::size_t block_size, bool symmetric=false)
Create a finite element from a Basix finite element.
Definition FiniteElement.cpp:97
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:287
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:574
const std::vector< std::vector< std::vector< int > > > & entity_closure_dofs() const noexcept
The local DOFs associated with the closure of each subentity of the cell.
Definition FiniteElement.cpp:244
std::span< const std::size_t > reference_value_shape() const noexcept
The reference value shape.
Definition FiniteElement.cpp:230
std::pair< std::vector< geometry_type >, std::array< std::size_t, 2 > > interpolation_operator() const
Definition FiniteElement.cpp:377
int block_size() const noexcept
Definition FiniteElement.cpp:263
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:639
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:418
bool interpolation_ident() const noexcept
Definition FiniteElement.cpp:351
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:622
bool operator==(const FiniteElement &e) const
Definition FiniteElement.cpp:199
bool map_ident() const noexcept
Definition FiniteElement.cpp:341
Finite element method functionality.
Definition assemble_matrix_impl.h:26
doftransform
DOF transformation type.
Definition FiniteElement.h:25
@ inverse_transpose
Transpose inverse.
CellType
Cell type identifier.
Definition cell_types.h:22