9 #include <basix/finite-element.h>
10 #include <dolfinx/common/types.h>
11 #include <dolfinx/mesh/cell_types.h>
16 #include <xtl/xspan.hpp>
18 struct ufc_finite_element;
29 explicit FiniteElement(
const ufc_finite_element& ufc_element);
82 std::
string family() const noexcept;
93 void tabulate(xt::xtensor<
double, 4>& values, const xt::xtensor<
double, 2>& X,
98 const xt::xtensor<
double, 3>& reference_values,
99 const xt::xtensor<
double, 3>& J,
100 const xtl::span<const
double>& detJ,
101 const xt::xtensor<
double, 3>& K) const;
112 std::
size_t hash() const noexcept;
151 template <typename T>
153 xtl::span<T> dofs)
const
157 throw std::runtime_error(
"No underlying element for interpolation. "
158 "Cannot interpolate mixed elements directly.");
161 _element->interpolate(tcb::make_span(dofs), tcb::make_span(values), _bs);
208 template <typename T>
209 std::function<
void(const xtl::span<T>&, const xtl::span<const std::uint32_t>&,
212 bool scalar_element = false)
const
217 return [](
const xtl::span<T>&,
const xtl::span<const std::uint32_t>&,
224 if (_sub_elements.size() != 0)
229 std::vector<std::function<void(
const xtl::span<T>&,
230 const xtl::span<const std::uint32_t>&,
232 sub_element_functions;
233 std::vector<int> dims;
234 for (std::size_t i = 0; i < _sub_elements.size(); ++i)
236 sub_element_functions.push_back(
237 _sub_elements[i]->get_dof_transformation_function<T>(inverse,
242 return [dims, sub_element_functions](
243 const xtl::span<T>& data,
244 const xtl::span<const std::uint32_t>& cell_info,
247 std::size_t start = 0;
248 for (std::size_t e = 0; e < sub_element_functions.size(); ++e)
250 const std::size_t width = dims[e] *
block_size;
251 sub_element_functions[e](data.subspan(start, width), cell_info,
257 else if (!scalar_element)
260 const std::function<void(
const xtl::span<T>&,
261 const xtl::span<const std::uint32_t>&,
263 sub_function = _sub_elements[0]->get_dof_transformation_function<T>(
267 [ebs, sub_function](
const xtl::span<T>& data,
268 const xtl::span<const std::uint32_t>& cell_info,
269 std::int32_t cell,
int data_block_size)
270 { sub_function(data, cell_info, cell, ebs * data_block_size); };
277 return [
this](
const xtl::span<T>& data,
278 const xtl::span<const std::uint32_t>& cell_info,
287 return [
this](
const xtl::span<T>& data,
288 const xtl::span<const std::uint32_t>& cell_info,
298 return [
this](
const xtl::span<T>& data,
299 const xtl::span<const std::uint32_t>& cell_info,
306 return [
this](
const xtl::span<T>& data,
307 const xtl::span<const std::uint32_t>& cell_info,
329 template <
typename T>
330 std::function<void(
const xtl::span<T>&,
const xtl::span<const std::uint32_t>&,
333 bool transpose =
false,
340 return [](
const xtl::span<T>&,
const xtl::span<const std::uint32_t>&,
346 else if (_sub_elements.size() != 0)
351 std::vector<std::function<void(
const xtl::span<T>&,
352 const xtl::span<const std::uint32_t>&,
354 sub_element_functions;
355 for (std::size_t i = 0; i < _sub_elements.size(); ++i)
357 sub_element_functions.push_back(
358 _sub_elements[i]->get_dof_transformation_to_transpose_function<T>(
359 inverse, transpose));
362 return [
this, sub_element_functions](
363 const xtl::span<T>& data,
364 const xtl::span<const std::uint32_t>& cell_info,
367 std::size_t start = 0;
368 for (std::size_t e = 0; e < sub_element_functions.size(); ++e)
370 const std::size_t width
371 = _sub_elements[e]->space_dimension() *
block_size;
372 sub_element_functions[e](data.subspan(start, width), cell_info,
378 else if (!scalar_element)
381 const std::function<void(
const xtl::span<T>&,
382 const xtl::span<const std::uint32_t>&,
384 sub_function = _sub_elements[0]->get_dof_transformation_function<T>(
387 sub_function](
const xtl::span<T>& data,
388 const xtl::span<const std::uint32_t>& cell_info,
389 std::int32_t cell,
int data_block_size)
392 const std::size_t dof_count = data.size() / data_block_size;
393 for (
int block = 0; block < data_block_size; ++block)
395 sub_function(data.subspan(block * dof_count, dof_count), cell_info,
406 return [
this](
const xtl::span<T>& data,
407 const xtl::span<const std::uint32_t>& cell_info,
416 return [
this](
const xtl::span<T>& data,
417 const xtl::span<const std::uint32_t>& cell_info,
429 return [
this](
const xtl::span<T>& data,
430 const xtl::span<const std::uint32_t>& cell_info,
439 return [
this](
const xtl::span<T>& data,
440 const xtl::span<const std::uint32_t>& cell_info,
454 template <
typename T>
456 std::uint32_t cell_permutation,
460 _element->apply_dof_transformation(data,
block_size, cell_permutation);
470 template <
typename T>
473 std::uint32_t cell_permutation,
477 _element->apply_inverse_transpose_dof_transformation(data,
block_size,
487 template <
typename T>
489 std::uint32_t cell_permutation,
493 _element->apply_transpose_dof_transformation(data,
block_size,
503 template <
typename T>
505 std::uint32_t cell_permutation,
509 _element->apply_inverse_dof_transformation(data,
block_size,
518 template <
typename T>
520 std::uint32_t cell_permutation,
524 _element->apply_dof_transformation_to_transpose(data,
block_size,
533 template <
typename T>
536 std::uint32_t cell_permutation,
540 _element->apply_inverse_dof_transformation_to_transpose(data,
block_size,
549 template <
typename T>
551 const xtl::span<T>& data, std::uint32_t cell_permutation,
555 _element->apply_transpose_dof_transformation_to_transpose(data,
block_size,
564 template <
typename T>
566 const xtl::span<T>& data, std::uint32_t cell_permutation,
570 _element->apply_inverse_transpose_dof_transformation_to_transpose(
592 template <
typename O,
typename P,
typename Q,
typename T,
typename S>
597 _element->map_pull_back_m(u, J, detJ, K, U);
604 void permute_dofs(
const xtl::span<std::int32_t>& doflist,
605 std::uint32_t cell_permutation)
const;
612 std::uint32_t cell_permutation)
const;
625 std::function<void(
const xtl::span<std::int32_t>&, std::uint32_t)>
627 bool scalar_element =
false)
const;
630 std::string _signature, _family;
634 int _tdim, _space_dim, _value_size, _reference_value_size;
637 std::vector<std::shared_ptr<const FiniteElement>> _sub_elements;
643 std::vector<int> _value_dimension;
650 bool _needs_dof_permutations;
651 bool _needs_dof_transformations;
654 std::unique_ptr<basix::FiniteElement> _element;
Finite Element, containing the dof layout on a reference element, and various methods for evaluating ...
Definition: FiniteElement.h:25
bool needs_dof_transformations() const noexcept
Check if DOF transformations are needed for this element.
Definition: FiniteElement.cpp:265
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.
Definition: FiniteElement.h:332
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...
Definition: FiniteElement.h:488
int space_dimension() const noexcept
Dimension of the finite element function space.
Definition: FiniteElement.cpp:179
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 transformat...
Definition: FiniteElement.h:472
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.
Definition: FiniteElement.cpp:288
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.
Definition: FiniteElement.h:211
FiniteElement & operator=(FiniteElement &&element)=default
Move assignment.
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.
Definition: FiniteElement.h:550
const std::vector< std::shared_ptr< const FiniteElement > > & sub_elements() const noexcept
Subelements (if any)
Definition: FiniteElement.cpp:228
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.
Definition: FiniteElement.cpp:183
std::string family() const noexcept
The finite element family.
Definition: FiniteElement.cpp:203
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.
Definition: FiniteElement.h:535
std::size_t hash() const noexcept
Return simple hash of the signature string.
Definition: FiniteElement.cpp:234
bool needs_dof_permutations() const noexcept
Check if DOF permutations are needed for this element.
Definition: FiniteElement.cpp:270
int value_dimension(int i) const
Return the dimension of the value space for axis i.
Definition: FiniteElement.cpp:195
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 ...
Definition: FiniteElement.cpp:253
void permute_dofs(const xtl::span< std::int32_t > &doflist, std::uint32_t cell_permutation) const
Permute the DOFs of the element.
Definition: FiniteElement.cpp:275
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.
Definition: FiniteElement.h:519
bool interpolation_ident() const noexcept
Check if interpolation into the finite element space is an identity operation given the evaluation on...
Definition: FiniteElement.cpp:247
std::string signature() const noexcept
String identifying the finite element.
Definition: FiniteElement.cpp:172
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.
Definition: FiniteElement.cpp:205
std::shared_ptr< const FiniteElement > extract_sub_element(const std::vector< int > &component) const
Extract sub finite element for component.
Definition: FiniteElement.cpp:237
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 t...
Definition: FiniteElement.h:504
FiniteElement(const ufc_finite_element &ufc_element)
Create finite element from UFC finite element.
Definition: FiniteElement.cpp:80
void unpermute_dofs(const xtl::span< std::int32_t > &doflist, std::uint32_t cell_permutation) const
Unpermute the DOFs of the element.
Definition: FiniteElement.cpp:281
int value_size() const noexcept
The value size, e.g. 1 for a scalar function, 2 for a 2D vector.
Definition: FiniteElement.cpp:181
int block_size() const noexcept
Block size of the finite element function space. For VectorElements and TensorElements,...
Definition: FiniteElement.cpp:193
int value_rank() const noexcept
Rank of the value space.
Definition: FiniteElement.cpp:188
FiniteElement & operator=(const FiniteElement &element)=delete
Copy assignment.
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 t...
Definition: FiniteElement.h:593
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.
Definition: FiniteElement.h:565
void apply_dof_transformation(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply DOF transformation to some data.
Definition: FiniteElement.h:455
mesh::CellType cell_shape() const noexcept
Cell shape.
Definition: FiniteElement.cpp:174
virtual ~FiniteElement()=default
Destructor.
int num_sub_elements() const noexcept
Get the number of sub elements (for a mixed element)
Definition: FiniteElement.cpp:222
constexpr void interpolate(const xt::xtensor< T, 2 > &values, xtl::span< T > dofs) const
Definition: FiniteElement.h:152
FiniteElement(FiniteElement &&element)=default
Move constructor.
FiniteElement(const FiniteElement &element)=delete
Copy constructor.
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.
Definition: FiniteElement.cpp:212
Finite element method functionality.
Definition: assemble_matrix_impl.h:23
CellType
Cell type identifier.
Definition: cell_types.h:22