8 #include "element-families.h"
12 #include "precompute.h"
13 #include "sobolev-spaces.h"
32 template <
typename T, std::
size_t d>
33 using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
34 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
35 template <
typename T, std::
size_t d>
37 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::mdarray<
38 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
43 std::array<std::vector<mdspan_t<const T, 2>>, 4>
44 to_mdspan(std::array<std::vector<mdarray_t<T, 2>>, 4>& x)
46 std::array<std::vector<mdspan_t<const T, 2>>, 4> x1;
47 for (std::size_t i = 0; i < x.size(); ++i)
48 for (std::size_t j = 0; j < x[i].size(); ++j)
49 x1[i].emplace_back(x[i][j].data(), x[i][j].extents());
57 std::array<std::vector<mdspan_t<const T, 4>>, 4>
58 to_mdspan(std::array<std::vector<mdarray_t<T, 4>>, 4>& M)
60 std::array<std::vector<mdspan_t<const T, 4>>, 4> M1;
61 for (std::size_t i = 0; i < M.size(); ++i)
62 for (std::size_t j = 0; j < M[i].size(); ++j)
63 M1[i].emplace_back(M[i][j].data(), M[i][j].extents());
71 std::array<std::vector<mdspan_t<const T, 2>>, 4>
72 to_mdspan(
const std::array<std::vector<std::vector<T>>, 4>& x,
73 const std::array<std::vector<std::array<std::size_t, 2>>, 4>& shape)
75 std::array<std::vector<mdspan_t<const T, 2>>, 4> x1;
76 for (std::size_t i = 0; i < x.size(); ++i)
77 for (std::size_t j = 0; j < x[i].size(); ++j)
78 x1[i].push_back(mdspan_t<const T, 2>(x[i][j].data(), shape[i][j]));
86 std::array<std::vector<mdspan_t<const T, 4>>, 4>
87 to_mdspan(
const std::array<std::vector<std::vector<T>>, 4>& M,
88 const std::array<std::vector<std::array<std::size_t, 4>>, 4>& shape)
90 std::array<std::vector<mdspan_t<const T, 4>>, 4> M1;
91 for (std::size_t i = 0; i < M.size(); ++i)
92 for (std::size_t j = 0; j < M[i].size(); ++j)
93 M1[i].push_back(mdspan_t<const T, 4>(M[i][j].data(), shape[i][j]));
103 template <
typename T, std::
size_t d>
121 template <std::
floating_po
int T>
122 std::tuple<std::array<std::vector<std::vector<T>>, 4>,
123 std::array<std::vector<std::array<std::size_t, 2>>, 4>,
124 std::array<std::vector<std::vector<T>>, 4>,
125 std::array<std::vector<std::array<std::size_t, 4>>, 4>>
128 std::size_t tdim, std::size_t value_size);
137 template <std::
floating_po
int F>
140 template <
typename T, std::
size_t d>
141 using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
142 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
321 const std::array<std::vector<mdspan_t<const F, 2>>, 4>&
x,
322 const std::array<std::vector<mdspan_t<const F, 4>>, 4>&
M,
352 std::size_t
hash()
const;
364 std::size_t num_points)
const
366 std::size_t ndsize = 1;
367 for (std::size_t i = 1; i <= nd; ++i)
368 ndsize *= (_cell_tdim + i);
369 for (std::size_t i = 1; i <= nd; ++i)
371 std::size_t vs = std::accumulate(_value_shape.begin(), _value_shape.end(),
372 1, std::multiplies{});
373 std::size_t ndofs = _coeffs.second[0];
374 return {ndsize, num_points, ndofs, vs};
398 std::pair<std::vector<F>, std::array<std::size_t, 4>>
399 tabulate(
int nd, impl::mdspan_t<const F, 2>
x)
const;
424 std::pair<std::vector<F>, std::array<std::size_t, 4>>
426 std::array<std::size_t, 2> shape)
const;
453 void tabulate(
int nd, impl::mdspan_t<const F, 2>
x,
454 mdspan_t<F, 4> basis)
const;
481 void tabulate(
int nd, std::span<const F>
x, std::array<std::size_t, 2> xshape,
482 std::span<F> basis)
const;
512 const std::vector<std::size_t>&
value_shape()
const {
return _value_shape; }
518 int dim()
const {
return _coeffs.second[0]; }
528 return _lagrange_variant;
553 return _dof_transformations_are_permutations;
559 return _dof_transformations_are_identity;
576 std::pair<std::vector<F>, std::array<std::size_t, 3>>
577 push_forward(impl::mdspan_t<const F, 3> U, impl::mdspan_t<const F, 3> J,
578 std::span<const F> detJ, impl::mdspan_t<const F, 3> K)
const;
587 std::pair<std::vector<F>, std::array<std::size_t, 3>>
588 pull_back(impl::mdspan_t<const F, 3> u, impl::mdspan_t<const F, 3> J,
589 std::span<const F> detJ, impl::mdspan_t<const F, 3> K)
const;
622 template <
typename O,
typename P,
typename Q,
typename R>
623 std::function<void(O&,
const P&,
const Q&, F,
const R&)>
map_fn()
const
627 case maps::type::identity:
628 return [](O& u,
const P& U,
const Q&, F,
const R&)
630 assert(U.extent(0) == u.extent(0));
631 assert(U.extent(1) == u.extent(1));
632 for (std::size_t i = 0; i < U.extent(0); ++i)
633 for (std::size_t j = 0; j < U.extent(1); ++j)
636 case maps::type::covariantPiola:
637 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
639 case maps::type::contravariantPiola:
640 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
642 case maps::type::doubleCovariantPiola:
643 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
645 case maps::type::doubleContravariantPiola:
646 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
649 throw std::runtime_error(
"Map not implemented");
660 const std::vector<std::vector<std::vector<int>>>&
entity_dofs()
const
676 return _e_closure_dofs;
760 std::pair<std::vector<F>, std::array<std::size_t, 3>>
767 std::map<cell::type, std::pair<std::vector<F>, std::array<std::size_t, 3>>>
770 return _entity_transformations;
795 void permute(std::span<std::int32_t> d, std::uint32_t cell_info)
const
797 if (!_dof_transformations_are_permutations)
799 throw std::runtime_error(
800 "The DOF transformations for this element are not permutations");
803 if (_dof_transformations_are_identity)
806 permute_data<std::int32_t, false>(d, 1, cell_info, _eperm);
826 void permute_inv(std::span<std::int32_t> d, std::uint32_t cell_info)
const
828 if (!_dof_transformations_are_permutations)
830 throw std::runtime_error(
831 "The DOF transformations for this element are not permutations");
834 if (_dof_transformations_are_identity)
837 permute_data<std::int32_t, true>(d, 1, cell_info, _eperm_inv);
856 std::uint32_t cell_info,
857 cell::type entity_type,
int entity_index)
const
861 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
863 std::uint32_t entity_info = 0;
869 entity_info = cell_info >> (face_start + entity_index) & 1;
872 entity_info = cell_info >> (3 * entity_index) & 7;
875 throw std::runtime_error(
"Unsupported cell dimension");
893 std::uint32_t cell_info,
895 int entity_index)
const
899 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
901 std::uint32_t entity_info;
907 entity_info = cell_info >> (face_start + entity_index) & 1;
910 entity_info = cell_info >> (3 * entity_index) & 7;
913 throw std::runtime_error(
"Unsupported cell dimension");
934 std::uint32_t entity_info,
937 if (!_dof_transformations_are_permutations)
939 throw std::runtime_error(
940 "The DOF transformations for this element are not permutations");
948 auto& perm = _subentity_closure_perm.at(entity_type);
956 else if (entity_dim == 2)
959 for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
972 throw std::runtime_error(
973 "Invalid dimension for permute_subentity_closure");
989 std::uint32_t entity_info,
992 if (!_dof_transformations_are_permutations)
994 throw std::runtime_error(
995 "The DOF transformations for this element are not permutations");
1000 if (entity_dim == 0)
1003 auto& perm = _subentity_closure_perm_inv.at(entity_type);
1004 if (entity_dim == 1)
1006 if (entity_info & 1)
1011 else if (entity_dim == 2)
1014 if (entity_info & 1)
1020 for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
1027 throw std::runtime_error(
1028 "Invalid dimension for permute_subentity_closure");
1062 template <
typename T>
1063 void T_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1075 template <
typename T>
1076 void Tt_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1089 template <
typename T>
1090 void Tt_inv_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1102 template <
typename T>
1103 void Tinv_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1115 template <
typename T>
1116 void Tt_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1127 template <
typename T>
1128 void T_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1140 template <
typename T>
1141 void Tinv_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1153 template <
typename T>
1162 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
points()
const
1218 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1229 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1232 return _dual_matrix;
1271 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
wcoeffs()
const
1281 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>, 4>&
1324 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>, 4>&
1333 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1351 return !_tensor_factors.empty();
1366 std::vector<std::vector<FiniteElement<F>>>
1370 throw std::runtime_error(
"Element has no tensor product representation.");
1371 return _tensor_factors;
1393 template <
typename T,
bool post>
1395 std::span<T> data,
int block_size, std::uint32_t cell_info,
1396 const std::map<
cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1399 using array2_t = std::pair<std::vector<F>, std::array<std::size_t, 2>>;
1400 using array3_t = std::pair<std::vector<F>, std::array<std::size_t, 3>>;
1402 = std::vector<std::pair<std::vector<std::size_t>, array2_t>>;
1410 template <
typename T,
bool post,
typename OP>
1412 transform_data(std::span<T> data,
int block_size, std::uint32_t cell_info,
1413 const std::map<cell::type, trans_data_t>& etrans, OP op)
const;
1422 std::size_t _cell_tdim;
1425 std::vector<std::vector<cell::type>> _cell_subentity_types;
1440 int _interpolation_nderivs;
1443 int _embedded_superdegree;
1446 int _embedded_subdegree;
1449 std::vector<std::size_t> _value_shape;
1462 std::pair<std::vector<F>, std::array<std::size_t, 2>> _coeffs;
1465 std::vector<std::vector<std::vector<int>>> _edofs;
1468 std::vector<std::vector<std::vector<int>>> _e_closure_dofs;
1471 std::map<cell::type, array3_t> _entity_transformations;
1478 std::pair<std::vector<F>, std::array<std::size_t, 2>> _points;
1482 std::array<std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>,
1487 std::pair<std::vector<F>, std::array<std::size_t, 2>> _matM;
1491 bool _dof_transformations_are_permutations;
1494 bool _dof_transformations_are_identity;
1499 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm;
1504 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm_inv;
1507 std::map<cell::type, trans_data_t> _etrans;
1510 std::map<cell::type, trans_data_t> _etransT;
1513 std::map<cell::type, trans_data_t> _etrans_inv;
1516 std::map<cell::type, trans_data_t> _etrans_invT;
1520 std::map<cell::type, std::vector<std::vector<std::size_t>>>
1521 _subentity_closure_perm;
1525 std::map<cell::type, std::vector<std::vector<std::size_t>>>
1526 _subentity_closure_perm_inv;
1530 bool _discontinuous;
1533 std::pair<std::vector<F>, std::array<std::size_t, 2>> _dual_matrix;
1541 std::vector<int> _dof_ordering;
1548 std::vector<std::vector<FiniteElement>> _tensor_factors;
1551 bool _interpolation_is_identity;
1555 std::pair<std::vector<F>, std::array<std::size_t, 2>> _wcoeffs;
1559 = std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>;
1560 std::array<array4_t, 4> _M;
1588 template <std::
floating_po
int T>
1590 cell::type cell_type,
const std::vector<std::size_t>& value_shape,
1591 impl::mdspan_t<const T, 2> wcoeffs,
1592 const std::array<std::vector<impl::mdspan_t<const T, 2>>, 4>& x,
1593 const std::array<std::vector<impl::mdspan_t<const T, 4>>, 4>& M,
1594 int interpolation_nderivs,
maps::type map_type,
1595 sobolev::space sobolev_space,
bool discontinuous,
int embedded_subdegree,
1609 template <std::
floating_po
int T>
1614 std::vector<int> dof_ordering = {});
1630 bool discontinuous);
1644 template <std::
floating_po
int T>
1645 std::vector<std::vector<FiniteElement<T>>>
1648 bool discontinuous, std::vector<int> dof_ordering);
1661 template <std::
floating_po
int T>
1672 template <std::
floating_po
int F>
1673 template <
typename T,
bool post>
1674 void FiniteElement<F>::permute_data(
1675 std::span<T> data,
int block_size, std::uint32_t cell_info,
1676 const std::map<
cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1679 if (_cell_tdim >= 2)
1683 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1687 auto& trans = eperm.at(cell::type::interval)[0];
1688 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1691 if (cell_info >> (face_start + e) & 1)
1699 if (_cell_tdim == 3)
1702 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1704 auto& trans = eperm.at(_cell_subentity_types[2][f]);
1707 if (!post and cell_info >> (3 * f) & 1)
1714 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1721 if (post and cell_info >> (3 * f) & 1)
1731 template <std::
floating_po
int F>
1732 template <
typename T,
bool post,
typename OP>
1733 void FiniteElement<F>::transform_data(
1734 std::span<T> data,
int block_size, std::uint32_t cell_info,
1735 const std::map<cell::type, trans_data_t>& etrans, OP op)
const
1737 if (_cell_tdim >= 2)
1741 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1743 for (
auto& edofs0 : _edofs[0])
1744 dofstart += edofs0.size();
1748 auto& [v_size_t, matrix] = etrans.at(cell::type::interval)[0];
1749 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1752 if (cell_info >> (face_start + e) & 1)
1754 op(std::span(v_size_t),
1755 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1756 dofstart, block_size);
1758 dofstart += _edofs[1][e].size();
1762 if (_cell_tdim == 3)
1765 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1767 auto& trans = etrans.at(_cell_subentity_types[2][f]);
1770 if (!post and cell_info >> (3 * f) & 1)
1772 const auto& m = trans[1];
1773 const auto& v_size_t = std::get<0>(m);
1774 const auto& matrix = std::get<1>(m);
1775 op(std::span(v_size_t),
1776 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1777 dofstart, block_size);
1781 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1783 const auto& m = trans[0];
1784 const auto& v_size_t = std::get<0>(m);
1785 const auto& matrix = std::get<1>(m);
1786 op(std::span(v_size_t),
1787 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1788 dofstart, block_size);
1792 if (post and cell_info >> (3 * f) & 1)
1794 const auto& m = trans[1];
1795 const auto& v_size_t = std::get<0>(m);
1796 const auto& matrix = std::get<1>(m);
1797 op(std::span(v_size_t),
1798 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1799 dofstart, block_size);
1802 dofstart += _edofs[2][f].size();
1808 template <std::
floating_po
int F>
1809 template <
typename T>
1811 std::uint32_t cell_info)
const
1813 if (_dof_transformations_are_identity)
1816 if (_dof_transformations_are_permutations)
1817 permute_data<T, false>(u, n, cell_info, _eperm);
1820 transform_data<T, false>(u, n, cell_info, _etrans,
1821 precompute::apply_matrix<F, T>);
1825 template <std::
floating_po
int F>
1826 template <
typename T>
1828 std::uint32_t cell_info)
const
1830 if (_dof_transformations_are_identity)
1832 else if (_dof_transformations_are_permutations)
1833 permute_data<T, true>(u, n, cell_info, _eperm_inv);
1836 transform_data<T, true>(u, n, cell_info, _etransT,
1837 precompute::apply_matrix<F, T>);
1841 template <std::
floating_po
int F>
1842 template <
typename T>
1844 std::uint32_t cell_info)
const
1846 if (_dof_transformations_are_identity)
1848 else if (_dof_transformations_are_permutations)
1849 permute_data<T, false>(u, n, cell_info, _eperm);
1852 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1853 precompute::apply_matrix<F, T>);
1857 template <std::
floating_po
int F>
1858 template <
typename T>
1860 std::uint32_t cell_info)
const
1862 if (_dof_transformations_are_identity)
1864 else if (_dof_transformations_are_permutations)
1865 permute_data<T, true>(u, n, cell_info, _eperm_inv);
1868 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1869 precompute::apply_matrix<F, T>);
1873 template <std::
floating_po
int F>
1874 template <
typename T>
1876 std::uint32_t cell_info)
const
1878 if (_dof_transformations_are_identity)
1880 else if (_dof_transformations_are_permutations)
1882 assert(u.size() % n == 0);
1883 const int step = u.size() / n;
1884 for (
int i = 0; i < n; ++i)
1886 std::span<T> dblock(u.data() + i * step, step);
1887 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1892 transform_data<T, false>(u, n, cell_info, _etrans,
1893 precompute::apply_tranpose_matrix_right<F, T>);
1897 template <std::
floating_po
int F>
1898 template <
typename T>
1900 std::uint32_t cell_info)
const
1902 if (_dof_transformations_are_identity)
1904 else if (_dof_transformations_are_permutations)
1906 assert(u.size() % n == 0);
1907 const int step = u.size() / n;
1908 for (
int i = 0; i < n; ++i)
1910 std::span<T> dblock(u.data() + i * step, step);
1911 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1916 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1917 precompute::apply_tranpose_matrix_right<F, T>);
1921 template <std::
floating_po
int F>
1922 template <
typename T>
1924 std::uint32_t cell_info)
const
1926 if (_dof_transformations_are_identity)
1928 else if (_dof_transformations_are_permutations)
1930 assert(u.size() % n == 0);
1931 const int step = u.size() / n;
1932 for (
int i = 0; i < n; ++i)
1934 std::span<T> dblock(u.data() + i * step, step);
1935 permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1940 transform_data<T, true>(u, n, cell_info, _etransT,
1941 precompute::apply_tranpose_matrix_right<F, T>);
1945 template <std::
floating_po
int F>
1946 template <
typename T>
1948 std::uint32_t cell_info)
const
1950 if (_dof_transformations_are_identity)
1952 else if (_dof_transformations_are_permutations)
1954 assert(u.size() % n == 0);
1955 const int step = u.size() / n;
1956 for (
int i = 0; i < n; ++i)
1958 std::span<T> dblock(u.data() + i * step, step);
1959 permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1964 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1965 precompute::apply_tranpose_matrix_right<F, T>);
A finite element.
Definition: finite-element.h:139
std::map< cell::type, std::pair< std::vector< F >, std::array< std::size_t, 3 > > > entity_transformations() const
Return the entity dof transformation matrices.
Definition: finite-element.h:768
void Tinv_apply_right(std::span< T > u, int n, std::uint32_t cell_info) const
Right(post)-apply the inverse of the operator applied by T_apply().
Definition: finite-element.h:1899
const std::vector< int > & dof_ordering() const
Get dof layout.
Definition: finite-element.h:1384
std::vector< std::vector< FiniteElement< F > > > get_tensor_product_representation() const
Get the tensor product representation of this element.
Definition: finite-element.h:1367
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & wcoeffs() const
Get the coefficients that define the polynomial set in terms of the orthonormal polynomials.
Definition: finite-element.h:1271
const std::array< std::vector< std::pair< std::vector< F >, std::array< std::size_t, 4 > > >, 4 > & M() const
Get the interpolation matrices for each subentity.
Definition: finite-element.h:1325
std::pair< std::vector< F >, std::array< std::size_t, 3 > > base_transformations() const
Get the base transformations.
Definition: finite-element.cpp:1632
void T_apply(std::span< T > u, int n, std::uint32_t cell_info) const
Transform basis functions from the reference element ordering and orientation to the globally consist...
Definition: finite-element.h:1810
void Tt_apply(std::span< T > u, int n, std::uint32_t cell_info) const
Apply the transpose of the operator applied by T_apply().
Definition: finite-element.h:1827
void permute_subentity_closure_inv(std::span< std::int32_t > d, std::uint32_t cell_info, cell::type entity_type, int entity_index) const
Perform the inverse of the operation applied by permute_subentity_closure().
Definition: finite-element.h:892
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & coefficient_matrix() const
Get the matrix of coefficients.
Definition: finite-element.h:1334
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & dual_matrix() const
Get the dual matrix.
Definition: finite-element.h:1230
bool dof_transformations_are_identity() const
Indicates is the dof transformations are all the identity.
Definition: finite-element.h:557
bool operator==(const FiniteElement &e) const
Check if two elements are the same.
Definition: finite-element.cpp:1451
void permute_subentity_closure_inv(std::span< std::int32_t > d, std::uint32_t entity_info, cell::type entity_type) const
Perform the inverse of the operation applied by permute_subentity_closure().
Definition: finite-element.h:988
int embedded_subdegree() const
Highest degree n such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this ...
Definition: finite-element.h:505
void T_apply_right(std::span< T > u, int n, std::uint32_t cell_info) const
Right(post)-apply the operator applied by T_apply().
Definition: finite-element.h:1923
FiniteElement(FiniteElement &&element)=default
Move constructor.
bool has_tensor_product_factorisation() const
Indicates whether or not this element can be represented as a product of elements defined on lower-di...
Definition: finite-element.h:1349
const std::vector< std::vector< std::vector< int > > > & entity_closure_dofs() const
Get the dofs on the closure of each topological entity: (vertices, edges, faces, cell) in that order.
Definition: finite-element.h:674
int interpolation_nderivs() const
The number of derivatives needed when interpolating.
Definition: finite-element.h:1381
std::pair< std::vector< F >, std::array< std::size_t, 4 > > tabulate(int nd, impl::mdspan_t< const F, 2 > x) const
Compute basis values and derivatives at set of points.
Definition: finite-element.cpp:1540
int embedded_superdegree() const
Lowest degree n such that the highest degree polynomial in this element is contained in a Lagrange (o...
Definition: finite-element.h:500
std::function< void(O &, const P &, const Q &, F, const R &)> map_fn() const
Return a function that performs the appropriate push-forward/pull-back for the element type.
Definition: finite-element.h:623
int dim() const
Dimension of the finite element space.
Definition: finite-element.h:518
FiniteElement(element::family family, cell::type cell_type, polyset::type poly_type, int degree, const std::vector< std::size_t > &value_shape, mdspan_t< const F, 2 > wcoeffs, const std::array< std::vector< mdspan_t< const F, 2 >>, 4 > &x, const std::array< std::vector< mdspan_t< const F, 4 >>, 4 > &M, int interpolation_nderivs, maps::type map_type, sobolev::space sobolev_space, bool discontinuous, int embedded_subdegree, int embedded_superdegree, element::lagrange_variant lvariant, element::dpc_variant dvariant, std::vector< int > dof_ordering={})
Construct a finite element.
void permute(std::span< std::int32_t > d, std::uint32_t cell_info) const
Permute indices associated with degree-of-freedoms on the reference element ordering to the globally ...
Definition: finite-element.h:795
std::size_t hash() const
Get a unique hash of this element.
Definition: finite-element.cpp:1490
void permute_inv(std::span< std::int32_t > d, std::uint32_t cell_info) const
Perform the inverse of the operation applied by permute().
Definition: finite-element.h:826
void Tt_apply_right(std::span< T > u, int n, std::uint32_t cell_info) const
Right(post)-apply the transpose of the operator applied by T_apply().
Definition: finite-element.h:1875
void Tt_inv_apply(std::span< T > u, int n, std::uint32_t cell_info) const
Apply the inverse transpose of the operator applied by T_apply().
Definition: finite-element.h:1843
void permute_subentity_closure(std::span< std::int32_t > d, std::uint32_t cell_info, cell::type entity_type, int entity_index) const
Permute indices associated with degree-of-freedoms on the closure of a sub-entity of the reference el...
Definition: finite-element.h:855
element::lagrange_variant lagrange_variant() const
Lagrange variant of the element.
Definition: finite-element.h:526
void Tt_inv_apply_right(std::span< T > u, int n, std::uint32_t cell_info) const
Right(post)-apply the transpose inverse of the operator applied by T_apply().
Definition: finite-element.h:1947
const std::vector< std::size_t > & value_shape() const
Element value tensor shape.
Definition: finite-element.h:512
int degree() const
Get the element polynomial degree.
Definition: finite-element.h:494
void permute_subentity_closure(std::span< std::int32_t > d, std::uint32_t entity_info, cell::type entity_type) const
Permute indices associated with degree-of-freedoms on the closure of a sub-entity of the reference el...
Definition: finite-element.h:933
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & points() const
Return the interpolation points.
Definition: finite-element.h:1162
void Tinv_apply(std::span< T > u, int n, std::uint32_t cell_info) const
Apply the inverse of the operator applied by T_apply().
Definition: finite-element.h:1859
std::array< std::size_t, 4 > tabulate_shape(std::size_t nd, std::size_t num_points) const
Array shape for tabulate basis values and derivatives at set of points.
Definition: finite-element.h:363
FiniteElement & operator=(const FiniteElement &element)=default
Assignment operator.
sobolev::space sobolev_space() const
Underlying Sobolev space for this element.
Definition: finite-element.h:541
FiniteElement(const FiniteElement &element)=default
Copy constructor.
polyset::type polyset_type() const
Get the element polyset type.
Definition: finite-element.h:490
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & interpolation_matrix() const
Return a matrix of weights interpolation.
Definition: finite-element.h:1219
std::pair< std::vector< F >, std::array< std::size_t, 3 > > pull_back(impl::mdspan_t< const F, 3 > u, impl::mdspan_t< const F, 3 > J, std::span< const F > detJ, impl::mdspan_t< const F, 3 > K) const
Map function values from a physical cell to the reference.
Definition: finite-element.cpp:1741
bool dof_transformations_are_permutations() const
Indicates if the degree-of-freedom transformations are all permutations.
Definition: finite-element.h:551
std::pair< std::vector< F >, std::array< std::size_t, 3 > > push_forward(impl::mdspan_t< const F, 3 > U, impl::mdspan_t< const F, 3 > J, std::span< const F > detJ, impl::mdspan_t< const F, 3 > K) const
Map function values from the reference to a physical cell.
Definition: finite-element.cpp:1701
const std::vector< std::vector< std::vector< int > > > & entity_dofs() const
Get the dofs on each topological entity: (vertices, edges, faces, cell) in that order.
Definition: finite-element.h:660
cell::type cell_type() const
Get the element cell type.
Definition: finite-element.h:486
bool interpolation_is_identity() const
Indicates whether or not the interpolation matrix for this element is an identity matrix.
Definition: finite-element.h:1378
const std::array< std::vector< std::pair< std::vector< F >, std::array< std::size_t, 2 > > >, 4 > & x() const
Get the interpolation points for each subentity.
Definition: finite-element.h:1282
FiniteElement & operator=(FiniteElement &&element)=default
Move assignment operator.
bool discontinuous() const
Indicates whether this element is the discontinuous variant.
Definition: finite-element.h:547
maps::type map_type() const
Map type for the element.
Definition: finite-element.h:537
element::dpc_variant dpc_variant() const
DPC variant of the element.
Definition: finite-element.h:533
element::family family() const
The finite element family.
Definition: finite-element.h:522
~FiniteElement()=default
Destructor.
type
Cell type.
Definition: cell.h:21
int topological_dimension(cell::type celltype)
Definition: cell.cpp:295
lagrange_variant
Variants of a Lagrange space that can be created.
Definition: element-families.h:12
impl::mdspan_t< T, d > mdspan_t
Typedef for mdspan.
Definition: finite-element.h:104
std::tuple< std::array< std::vector< std::vector< T > >, 4 >, std::array< std::vector< std::array< std::size_t, 2 > >, 4 >, std::array< std::vector< std::vector< T > >, 4 >, std::array< std::vector< std::array< std::size_t, 4 > >, 4 > > make_discontinuous(const std::array< std::vector< mdspan_t< const T, 2 >>, 4 > &x, const std::array< std::vector< mdspan_t< const T, 4 >>, 4 > &M, std::size_t tdim, std::size_t value_size)
Definition: finite-element.cpp:521
dpc_variant
Definition: element-families.h:32
family
Available element families.
Definition: element-families.h:45
void covariant_piola(O &&r, const P &U, const Q &, double, const R &K)
Covariant Piola map.
Definition: maps.h:61
void contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Contravariant Piola map.
Definition: maps.h:81
void double_contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Double contravariant Piola map.
Definition: maps.h:135
void double_covariant_piola(O &&r, const P &U, const Q &J, double, const R &K)
Double covariant Piola map.
Definition: maps.h:103
type
Map type.
Definition: maps.h:39
type
Cell type.
Definition: polyset.h:136
void apply_permutation(std::span< const std::size_t > perm, std::span< E > data, std::size_t offset=0, std::size_t n=1)
Apply a (precomputed) permutation .
Definition: precompute.h:144
void apply_permutation_mapped(std::span< const std::size_t > perm, std::span< E > data, std::span< const int > emap, std::size_t n=1)
Permutation of mapped data.
Definition: precompute.h:154
space
Sobolev space type.
Definition: sobolev-spaces.h:13
Basix: FEniCS runtime basis evaluation library.
Definition: cell.h:17
FiniteElement< T > create_element(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous, std::vector< int > dof_ordering={})
Definition: finite-element.cpp:194
std::vector< int > tp_dof_ordering(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous)
Definition: finite-element.cpp:400
std::vector< std::vector< FiniteElement< T > > > tp_factors(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous, std::vector< int > dof_ordering)
Definition: finite-element.cpp:350
std::string version()
Definition: finite-element.cpp:1779
FiniteElement< T > create_tp_element(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous)
Definition: finite-element.cpp:331
FiniteElement< T > create_custom_element(cell::type cell_type, const std::vector< std::size_t > &value_shape, impl::mdspan_t< const T, 2 > wcoeffs, const std::array< std::vector< impl::mdspan_t< const T, 2 >>, 4 > &x, const std::array< std::vector< impl::mdspan_t< const T, 4 >>, 4 > &M, int interpolation_nderivs, maps::type map_type, sobolev::space sobolev_space, bool discontinuous, int embedded_subdegree, int embedded_superdegree, polyset::type poly_type)
Definition: finite-element.cpp:610