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 = md::mdspan<T, md::dextents<std::size_t, d>>;
34 template <
typename T, std::
size_t d>
36 = md::MDSPAN_IMPL_PROPOSED_NAMESPACE::mdarray<T,
37 md::dextents<std::size_t, d>>;
42 std::array<std::vector<mdspan_t<const T, 2>>, 4>
43 to_mdspan(std::array<std::vector<mdarray_t<T, 2>>, 4>& x)
45 std::array<std::vector<mdspan_t<const T, 2>>, 4> x1;
46 for (std::size_t i = 0; i < x.size(); ++i)
47 for (std::size_t j = 0; j < x[i].size(); ++j)
48 x1[i].emplace_back(x[i][j].data(), x[i][j].extents());
56 std::array<std::vector<mdspan_t<const T, 4>>, 4>
57 to_mdspan(std::array<std::vector<mdarray_t<T, 4>>, 4>& M)
59 std::array<std::vector<mdspan_t<const T, 4>>, 4> M1;
60 for (std::size_t i = 0; i < M.size(); ++i)
61 for (std::size_t j = 0; j < M[i].size(); ++j)
62 M1[i].emplace_back(M[i][j].data(), M[i][j].extents());
70 std::array<std::vector<mdspan_t<const T, 2>>, 4>
71 to_mdspan(
const std::array<std::vector<std::vector<T>>, 4>& x,
72 const std::array<std::vector<std::array<std::size_t, 2>>, 4>& shape)
74 std::array<std::vector<mdspan_t<const T, 2>>, 4> x1;
75 for (std::size_t i = 0; i < x.size(); ++i)
76 for (std::size_t j = 0; j < x[i].size(); ++j)
77 x1[i].push_back(mdspan_t<const T, 2>(x[i][j].data(), shape[i][j]));
85 std::array<std::vector<mdspan_t<const T, 4>>, 4>
86 to_mdspan(
const std::array<std::vector<std::vector<T>>, 4>& M,
87 const std::array<std::vector<std::array<std::size_t, 4>>, 4>& shape)
89 std::array<std::vector<mdspan_t<const T, 4>>, 4> M1;
90 for (std::size_t i = 0; i < M.size(); ++i)
91 for (std::size_t j = 0; j < M[i].size(); ++j)
92 M1[i].push_back(mdspan_t<const T, 4>(M[i][j].data(), shape[i][j]));
102 template <
typename T, std::
size_t d>
120 template <std::
floating_po
int T>
121 std::tuple<std::array<std::vector<std::vector<T>>, 4>,
122 std::array<std::vector<std::array<std::size_t, 2>>, 4>,
123 std::array<std::vector<std::vector<T>>, 4>,
124 std::array<std::vector<std::array<std::size_t, 4>>, 4>>
127 std::size_t tdim, std::size_t value_size);
136 template <std::
floating_po
int F>
139 template <
typename T, std::
size_t d>
140 using mdspan_t = md::mdspan<T, md::dextents<std::size_t, d>>;
322 const std::array<std::vector<mdspan_t<const F, 2>>, 4>&
x,
323 const std::array<std::vector<mdspan_t<const F, 4>>, 4>&
M,
353 std::size_t
hash()
const;
365 std::size_t num_points)
const
367 std::size_t ndsize = 1;
368 for (std::size_t i = 1; i <= nd; ++i)
369 ndsize *= (_cell_tdim + i);
370 for (std::size_t i = 1; i <= nd; ++i)
372 std::size_t vs = std::accumulate(_value_shape.begin(), _value_shape.end(),
373 1, std::multiplies{});
374 std::size_t ndofs = _coeffs.second[0];
375 return {ndsize, num_points, ndofs, vs};
399 std::pair<std::vector<F>, std::array<std::size_t, 4>>
400 tabulate(
int nd, impl::mdspan_t<const F, 2>
x)
const;
425 std::pair<std::vector<F>, std::array<std::size_t, 4>>
427 std::array<std::size_t, 2> shape)
const;
454 void tabulate(
int nd, impl::mdspan_t<const F, 2>
x,
455 mdspan_t<F, 4> basis)
const;
482 void tabulate(
int nd, std::span<const F>
x, std::array<std::size_t, 2> xshape,
483 std::span<F> basis)
const;
513 const std::vector<std::size_t>&
value_shape()
const {
return _value_shape; }
519 int dim()
const {
return _coeffs.second[0]; }
529 return _lagrange_variant;
554 return _dof_transformations_are_permutations;
560 return _dof_transformations_are_identity;
577 std::pair<std::vector<F>, std::array<std::size_t, 3>>
578 push_forward(impl::mdspan_t<const F, 3> U, impl::mdspan_t<const F, 3> J,
579 std::span<const F> detJ, impl::mdspan_t<const F, 3> K)
const;
588 std::pair<std::vector<F>, std::array<std::size_t, 3>>
589 pull_back(impl::mdspan_t<const F, 3> u, impl::mdspan_t<const F, 3> J,
590 std::span<const F> detJ, impl::mdspan_t<const F, 3> K)
const;
623 template <
typename O,
typename P,
typename Q,
typename R>
624 std::function<void(O&,
const P&,
const Q&, F,
const R&)>
map_fn()
const
628 case maps::type::identity:
629 return [](O& u,
const P& U,
const Q&, F,
const R&)
631 assert(U.extent(0) == u.extent(0));
632 assert(U.extent(1) == u.extent(1));
633 for (std::size_t i = 0; i < U.extent(0); ++i)
634 for (std::size_t j = 0; j < U.extent(1); ++j)
637 case maps::type::covariantPiola:
638 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
640 case maps::type::contravariantPiola:
641 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
643 case maps::type::doubleCovariantPiola:
644 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
646 case maps::type::doubleContravariantPiola:
647 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
650 throw std::runtime_error(
"Map not implemented");
661 const std::vector<std::vector<std::vector<int>>>&
entity_dofs()
const
677 return _e_closure_dofs;
761 std::pair<std::vector<F>, std::array<std::size_t, 3>>
768 std::map<cell::type, std::pair<std::vector<F>, std::array<std::size_t, 3>>>
771 return _entity_transformations;
796 void permute(std::span<std::int32_t> d, std::uint32_t cell_info)
const
798 if (!_dof_transformations_are_permutations)
800 throw std::runtime_error(
801 "The DOF transformations for this element are not permutations");
804 if (_dof_transformations_are_identity)
807 permute_data<std::int32_t, false>(d, 1, cell_info, _eperm);
827 void permute_inv(std::span<std::int32_t> d, std::uint32_t cell_info)
const
829 if (!_dof_transformations_are_permutations)
831 throw std::runtime_error(
832 "The DOF transformations for this element are not permutations");
835 if (_dof_transformations_are_identity)
838 permute_data<std::int32_t, true>(d, 1, cell_info, _eperm_inv);
857 std::uint32_t cell_info,
858 cell::type entity_type,
int entity_index)
const
862 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
864 std::uint32_t entity_info = 0;
870 entity_info = cell_info >> (face_start + entity_index) & 1;
873 entity_info = cell_info >> (3 * entity_index) & 7;
876 throw std::runtime_error(
"Unsupported cell dimension");
894 std::uint32_t cell_info,
896 int entity_index)
const
900 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
902 std::uint32_t entity_info;
908 entity_info = cell_info >> (face_start + entity_index) & 1;
911 entity_info = cell_info >> (3 * entity_index) & 7;
914 throw std::runtime_error(
"Unsupported cell dimension");
935 std::uint32_t entity_info,
938 if (!_dof_transformations_are_permutations)
940 throw std::runtime_error(
941 "The DOF transformations for this element are not permutations");
949 auto& perm = _subentity_closure_perm.at(entity_type);
957 else if (entity_dim == 2)
960 for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
973 throw std::runtime_error(
974 "Invalid dimension for permute_subentity_closure");
990 std::uint32_t entity_info,
993 if (!_dof_transformations_are_permutations)
995 throw std::runtime_error(
996 "The DOF transformations for this element are not permutations");
1001 if (entity_dim == 0)
1004 auto& perm = _subentity_closure_perm_inv.at(entity_type);
1005 if (entity_dim == 1)
1007 if (entity_info & 1)
1012 else if (entity_dim == 2)
1015 if (entity_info & 1)
1021 for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
1028 throw std::runtime_error(
1029 "Invalid dimension for permute_subentity_closure");
1063 template <
typename T>
1064 void T_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1076 template <
typename T>
1077 void Tt_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1090 template <
typename T>
1091 void Tt_inv_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1103 template <
typename T>
1104 void Tinv_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1116 template <
typename T>
1117 void Tt_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1128 template <
typename T>
1129 void T_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1141 template <
typename T>
1142 void Tinv_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1154 template <
typename T>
1163 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
points()
const
1219 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1230 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1233 return _dual_matrix;
1272 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
wcoeffs()
const
1282 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>, 4>&
1325 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>, 4>&
1334 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1352 return !_tensor_factors.empty();
1367 std::vector<std::vector<FiniteElement<F>>>
1371 throw std::runtime_error(
"Element has no tensor product representation.");
1372 return _tensor_factors;
1394 template <
typename T,
bool post>
1396 std::span<T> data,
int block_size, std::uint32_t cell_info,
1397 const std::map<
cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1400 using array2_t = std::pair<std::vector<F>, std::array<std::size_t, 2>>;
1401 using array3_t = std::pair<std::vector<F>, std::array<std::size_t, 3>>;
1403 = std::vector<std::pair<std::vector<std::size_t>, array2_t>>;
1411 template <
typename T,
bool post,
typename OP>
1413 transform_data(std::span<T> data,
int block_size, std::uint32_t cell_info,
1414 const std::map<cell::type, trans_data_t>& etrans, OP op)
const;
1423 std::size_t _cell_tdim;
1426 std::vector<std::vector<cell::type>> _cell_subentity_types;
1441 int _interpolation_nderivs;
1444 int _embedded_superdegree;
1447 int _embedded_subdegree;
1450 std::vector<std::size_t> _value_shape;
1463 std::pair<std::vector<F>, std::array<std::size_t, 2>> _coeffs;
1466 std::vector<std::vector<std::vector<int>>> _edofs;
1469 std::vector<std::vector<std::vector<int>>> _e_closure_dofs;
1472 std::map<cell::type, array3_t> _entity_transformations;
1479 std::pair<std::vector<F>, std::array<std::size_t, 2>> _points;
1483 std::array<std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>,
1488 std::pair<std::vector<F>, std::array<std::size_t, 2>> _matM;
1492 bool _dof_transformations_are_permutations;
1495 bool _dof_transformations_are_identity;
1500 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm;
1505 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm_inv;
1508 std::map<cell::type, trans_data_t> _etrans;
1511 std::map<cell::type, trans_data_t> _etransT;
1514 std::map<cell::type, trans_data_t> _etrans_inv;
1517 std::map<cell::type, trans_data_t> _etrans_invT;
1521 std::map<cell::type, std::vector<std::vector<std::size_t>>>
1522 _subentity_closure_perm;
1526 std::map<cell::type, std::vector<std::vector<std::size_t>>>
1527 _subentity_closure_perm_inv;
1531 bool _discontinuous;
1534 std::pair<std::vector<F>, std::array<std::size_t, 2>> _dual_matrix;
1542 std::vector<int> _dof_ordering;
1549 std::vector<std::vector<FiniteElement>> _tensor_factors;
1552 bool _interpolation_is_identity;
1556 std::pair<std::vector<F>, std::array<std::size_t, 2>> _wcoeffs;
1560 = std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>;
1561 std::array<array4_t, 4> _M;
1589 template <std::
floating_po
int T>
1591 cell::type cell_type,
const std::vector<std::size_t>& value_shape,
1592 impl::mdspan_t<const T, 2> wcoeffs,
1593 const std::array<std::vector<impl::mdspan_t<const T, 2>>, 4>& x,
1594 const std::array<std::vector<impl::mdspan_t<const T, 4>>, 4>& M,
1595 int interpolation_nderivs,
maps::type map_type,
1596 sobolev::space sobolev_space,
bool discontinuous,
int embedded_subdegree,
1610 template <std::
floating_po
int T>
1615 std::vector<int> dof_ordering = {});
1629 std::optional<std::vector<int>>
1647 bool discontinuous);
1662 template <std::
floating_po
int T>
1663 std::optional<std::vector<std::vector<FiniteElement<T>>>>
1666 bool discontinuous,
const std::vector<int>& dof_ordering);
1679 template <std::
floating_po
int T>
1690 template <std::
floating_po
int F>
1691 template <
typename T,
bool post>
1692 void FiniteElement<F>::permute_data(
1693 std::span<T> data,
int block_size, std::uint32_t cell_info,
1694 const std::map<
cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1697 if (_cell_tdim >= 2)
1701 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1705 auto& trans = eperm.at(cell::type::interval)[0];
1706 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1709 if (cell_info >> (face_start + e) & 1)
1717 if (_cell_tdim == 3)
1720 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1722 auto& trans = eperm.at(_cell_subentity_types[2][f]);
1725 if (!post and cell_info >> (3 * f) & 1)
1732 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1739 if (post and cell_info >> (3 * f) & 1)
1749 template <std::
floating_po
int F>
1750 template <
typename T,
bool post,
typename OP>
1751 void FiniteElement<F>::transform_data(
1752 std::span<T> data,
int block_size, std::uint32_t cell_info,
1753 const std::map<cell::type, trans_data_t>& etrans, OP op)
const
1755 if (_cell_tdim >= 2)
1759 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1761 for (
auto& edofs0 : _edofs[0])
1762 dofstart += edofs0.size();
1766 auto& [v_size_t, matrix] = etrans.at(cell::type::interval)[0];
1767 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1770 if (cell_info >> (face_start + e) & 1)
1772 op(std::span(v_size_t),
1773 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1774 dofstart, block_size);
1776 dofstart += _edofs[1][e].size();
1780 if (_cell_tdim == 3)
1783 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1785 auto& trans = etrans.at(_cell_subentity_types[2][f]);
1788 if (!post and cell_info >> (3 * f) & 1)
1790 const auto& m = trans[1];
1791 const auto& v_size_t = std::get<0>(m);
1792 const auto& matrix = std::get<1>(m);
1793 op(std::span(v_size_t),
1794 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1795 dofstart, block_size);
1799 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1801 const auto& m = trans[0];
1802 const auto& v_size_t = std::get<0>(m);
1803 const auto& matrix = std::get<1>(m);
1804 op(std::span(v_size_t),
1805 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1806 dofstart, block_size);
1810 if (post and cell_info >> (3 * f) & 1)
1812 const auto& m = trans[1];
1813 const auto& v_size_t = std::get<0>(m);
1814 const auto& matrix = std::get<1>(m);
1815 op(std::span(v_size_t),
1816 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1817 dofstart, block_size);
1820 dofstart += _edofs[2][f].size();
1826 template <std::
floating_po
int F>
1827 template <
typename T>
1829 std::uint32_t cell_info)
const
1831 if (_dof_transformations_are_identity)
1834 if (_dof_transformations_are_permutations)
1835 permute_data<T, false>(u, n, cell_info, _eperm);
1838 transform_data<T, false>(u, n, cell_info, _etrans,
1839 precompute::apply_matrix<F, T>);
1843 template <std::
floating_po
int F>
1844 template <
typename T>
1846 std::uint32_t cell_info)
const
1848 if (_dof_transformations_are_identity)
1850 else if (_dof_transformations_are_permutations)
1851 permute_data<T, true>(u, n, cell_info, _eperm_inv);
1854 transform_data<T, true>(u, n, cell_info, _etransT,
1855 precompute::apply_matrix<F, T>);
1859 template <std::
floating_po
int F>
1860 template <
typename T>
1862 std::uint32_t cell_info)
const
1864 if (_dof_transformations_are_identity)
1866 else if (_dof_transformations_are_permutations)
1867 permute_data<T, false>(u, n, cell_info, _eperm);
1870 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1871 precompute::apply_matrix<F, T>);
1875 template <std::
floating_po
int F>
1876 template <
typename T>
1878 std::uint32_t cell_info)
const
1880 if (_dof_transformations_are_identity)
1882 else if (_dof_transformations_are_permutations)
1883 permute_data<T, true>(u, n, cell_info, _eperm_inv);
1886 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1887 precompute::apply_matrix<F, T>);
1891 template <std::
floating_po
int F>
1892 template <
typename T>
1894 std::uint32_t cell_info)
const
1896 if (_dof_transformations_are_identity)
1898 else if (_dof_transformations_are_permutations)
1900 assert(u.size() % n == 0);
1901 const int step = u.size() / n;
1902 for (
int i = 0; i < n; ++i)
1904 std::span<T> dblock(u.data() + i * step, step);
1905 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1910 transform_data<T, false>(u, n, cell_info, _etrans,
1911 precompute::apply_tranpose_matrix_right<F, T>);
1915 template <std::
floating_po
int F>
1916 template <
typename T>
1918 std::uint32_t cell_info)
const
1920 if (_dof_transformations_are_identity)
1922 else if (_dof_transformations_are_permutations)
1924 assert(u.size() % n == 0);
1925 const int step = u.size() / n;
1926 for (
int i = 0; i < n; ++i)
1928 std::span<T> dblock(u.data() + i * step, step);
1929 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1934 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1935 precompute::apply_tranpose_matrix_right<F, T>);
1939 template <std::
floating_po
int F>
1940 template <
typename T>
1942 std::uint32_t cell_info)
const
1944 if (_dof_transformations_are_identity)
1946 else if (_dof_transformations_are_permutations)
1948 assert(u.size() % n == 0);
1949 const int step = u.size() / n;
1950 for (
int i = 0; i < n; ++i)
1952 std::span<T> dblock(u.data() + i * step, step);
1953 permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1958 transform_data<T, true>(u, n, cell_info, _etransT,
1959 precompute::apply_tranpose_matrix_right<F, T>);
1963 template <std::
floating_po
int F>
1964 template <
typename T>
1966 std::uint32_t cell_info)
const
1968 if (_dof_transformations_are_identity)
1970 else if (_dof_transformations_are_permutations)
1972 assert(u.size() % n == 0);
1973 const int step = u.size() / n;
1974 for (
int i = 0; i < n; ++i)
1976 std::span<T> dblock(u.data() + i * step, step);
1977 permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1982 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1983 precompute::apply_tranpose_matrix_right<F, T>);
A finite element.
Definition: finite-element.h:138
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:769
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:1917
const std::vector< int > & dof_ordering() const
Get dof layout.
Definition: finite-element.h:1385
std::vector< std::vector< FiniteElement< F > > > get_tensor_product_representation() const
Get the tensor product representation of this element.
Definition: finite-element.h:1368
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:1272
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:1326
std::pair< std::vector< F >, std::array< std::size_t, 3 > > base_transformations() const
Get the base transformations.
Definition: finite-element.cpp:1881
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:1828
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:1845
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:893
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & coefficient_matrix() const
Get the matrix of coefficients.
Definition: finite-element.h:1335
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & dual_matrix() const
Get the dual matrix.
Definition: finite-element.h:1231
bool dof_transformations_are_identity() const
Indicates is the dof transformations are all the identity.
Definition: finite-element.h:558
bool operator==(const FiniteElement &e) const
Check if two elements are the same.
Definition: finite-element.cpp:1703
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:989
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:506
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:1941
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:1350
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:675
int interpolation_nderivs() const
The number of derivatives needed when interpolating.
Definition: finite-element.h:1382
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:1789
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:501
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:624
int dim() const
Dimension of the finite element space.
Definition: finite-element.h:519
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:796
std::size_t hash() const
Get a unique hash of this element.
Definition: finite-element.cpp:1742
F scalar_type
Scalar type.
Definition: finite-element.h:144
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:827
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:1893
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:1861
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:856
element::lagrange_variant lagrange_variant() const
Lagrange variant of the element.
Definition: finite-element.h:527
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:1965
const std::vector< std::size_t > & value_shape() const
Element value tensor shape.
Definition: finite-element.h:513
int degree() const
Get the element polynomial degree.
Definition: finite-element.h:495
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:934
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & points() const
Return the interpolation points.
Definition: finite-element.h:1163
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:1877
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:364
FiniteElement & operator=(const FiniteElement &element)=default
Assignment operator.
sobolev::space sobolev_space() const
Underlying Sobolev space for this element.
Definition: finite-element.h:542
FiniteElement(const FiniteElement &element)=default
Copy constructor.
polyset::type polyset_type() const
Get the element polyset type.
Definition: finite-element.h:491
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:1220
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:1986
bool dof_transformations_are_permutations() const
Indicates if the degree-of-freedom transformations are all permutations.
Definition: finite-element.h:552
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:1950
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:661
cell::type cell_type() const
Get the element cell type.
Definition: finite-element.h:487
bool interpolation_is_identity() const
Indicates whether or not the interpolation matrix for this element is an identity matrix.
Definition: finite-element.h:1379
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:1283
FiniteElement & operator=(FiniteElement &&element)=default
Move assignment operator.
bool discontinuous() const
Indicates whether this element is the discontinuous variant.
Definition: finite-element.h:548
maps::type map_type() const
Map type for the element.
Definition: finite-element.h:538
element::dpc_variant dpc_variant() const
DPC variant of the element.
Definition: finite-element.h:534
element::family family() const
The finite element family.
Definition: finite-element.h:523
~FiniteElement()=default
Destructor.
type
Cell type.
Definition: cell.h:21
int topological_dimension(cell::type celltype)
Definition: cell.cpp:294
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:103
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:816
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:62
void contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Contravariant Piola map.
Definition: maps.h:82
void double_contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Double contravariant Piola map.
Definition: maps.h:132
void double_covariant_piola(O &&r, const P &U, const Q &, double, const R &K)
Double covariant Piola map.
Definition: maps.h:104
type
Map type.
Definition: maps.h:40
type
Cell type.
Definition: polyset.h:137
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:137
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:147
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:192
std::optional< 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, const std::vector< int > &dof_ordering)
Definition: finite-element.cpp:353
std::string version()
Definition: finite-element.cpp:2020
std::vector< int > lex_dof_ordering(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous)
Definition: finite-element.cpp:509
std::optional< 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:398
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:329
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:901