8 #include "element-families.h"
12 #include "precompute.h"
13 #include "sobolev-spaces.h"
31 template <
typename T, std::
size_t d>
32 using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
33 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
34 template <
typename T, std::
size_t d>
36 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::mdarray<
37 T, MDSPAN_IMPL_STANDARD_NAMESPACE::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 = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
141 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
323 const std::array<std::vector<mdspan_t<const F, 2>>, 4>&
x,
324 const std::array<std::vector<mdspan_t<const F, 4>>, 4>&
M,
354 std::size_t
hash()
const;
366 std::size_t num_points)
const
368 std::size_t ndsize = 1;
369 for (std::size_t i = 1; i <= nd; ++i)
370 ndsize *= (_cell_tdim + i);
371 for (std::size_t i = 1; i <= nd; ++i)
373 std::size_t vs = std::accumulate(_value_shape.begin(), _value_shape.end(),
374 1, std::multiplies{});
375 std::size_t ndofs = _coeffs.second[0];
376 return {ndsize, num_points, ndofs, vs};
400 std::pair<std::vector<F>, std::array<std::size_t, 4>>
401 tabulate(
int nd, impl::mdspan_t<const F, 2>
x)
const;
426 std::pair<std::vector<F>, std::array<std::size_t, 4>>
428 std::array<std::size_t, 2> shape)
const;
455 void tabulate(
int nd, impl::mdspan_t<const F, 2>
x,
456 mdspan_t<F, 4> basis)
const;
483 void tabulate(
int nd, std::span<const F>
x, std::array<std::size_t, 2> xshape,
484 std::span<F> basis)
const;
514 const std::vector<std::size_t>&
value_shape()
const {
return _value_shape; }
520 int dim()
const {
return _coeffs.second[0]; }
530 return _lagrange_variant;
555 return _dof_transformations_are_permutations;
561 return _dof_transformations_are_identity;
578 std::pair<std::vector<F>, std::array<std::size_t, 3>>
579 push_forward(impl::mdspan_t<const F, 3> U, impl::mdspan_t<const F, 3> J,
580 std::span<const F> detJ, impl::mdspan_t<const F, 3> K)
const;
589 std::pair<std::vector<F>, std::array<std::size_t, 3>>
590 pull_back(impl::mdspan_t<const F, 3> u, impl::mdspan_t<const F, 3> J,
591 std::span<const F> detJ, impl::mdspan_t<const F, 3> K)
const;
624 template <
typename O,
typename P,
typename Q,
typename R>
625 std::function<void(O&,
const P&,
const Q&, F,
const R&)>
map_fn()
const
629 case maps::type::identity:
630 return [](O& u,
const P& U,
const Q&, F,
const R&)
632 assert(U.extent(0) == u.extent(0));
633 assert(U.extent(1) == u.extent(1));
634 for (std::size_t i = 0; i < U.extent(0); ++i)
635 for (std::size_t j = 0; j < U.extent(1); ++j)
638 case maps::type::covariantPiola:
639 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
641 case maps::type::contravariantPiola:
642 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
644 case maps::type::doubleCovariantPiola:
645 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
647 case maps::type::doubleContravariantPiola:
648 return [](O& u,
const P& U,
const Q& J, F detJ,
const R& K)
651 throw std::runtime_error(
"Map not implemented");
662 const std::vector<std::vector<std::vector<int>>>&
entity_dofs()
const
678 return _e_closure_dofs;
762 std::pair<std::vector<F>, std::array<std::size_t, 3>>
769 std::map<cell::type, std::pair<std::vector<F>, std::array<std::size_t, 3>>>
772 return _entity_transformations;
797 void permute(std::span<std::int32_t> d, std::uint32_t cell_info)
const
799 if (!_dof_transformations_are_permutations)
801 throw std::runtime_error(
802 "The DOF transformations for this element are not permutations");
805 if (_dof_transformations_are_identity)
808 permute_data<std::int32_t, false>(d, 1, cell_info, _eperm);
828 void permute_inv(std::span<std::int32_t> d, std::uint32_t cell_info)
const
830 if (!_dof_transformations_are_permutations)
832 throw std::runtime_error(
833 "The DOF transformations for this element are not permutations");
836 if (_dof_transformations_are_identity)
839 permute_data<std::int32_t, true>(d, 1, cell_info, _eperm_inv);
858 std::uint32_t cell_info,
859 cell::type entity_type,
int entity_index)
const
863 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
865 std::uint32_t entity_info = 0;
871 entity_info = cell_info >> (face_start + entity_index) & 1;
874 entity_info = cell_info >> (3 * entity_index) & 7;
877 throw std::runtime_error(
"Unsupported cell dimension");
895 std::uint32_t cell_info,
897 int entity_index)
const
901 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
903 std::uint32_t entity_info;
909 entity_info = cell_info >> (face_start + entity_index) & 1;
912 entity_info = cell_info >> (3 * entity_index) & 7;
915 throw std::runtime_error(
"Unsupported cell dimension");
936 std::uint32_t entity_info,
939 if (!_dof_transformations_are_permutations)
941 throw std::runtime_error(
942 "The DOF transformations for this element are not permutations");
950 auto& perm = _subentity_closure_perm.at(entity_type);
958 else if (entity_dim == 2)
961 for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
974 throw std::runtime_error(
975 "Invalid dimension for permute_subentity_closure");
991 std::uint32_t entity_info,
994 if (!_dof_transformations_are_permutations)
996 throw std::runtime_error(
997 "The DOF transformations for this element are not permutations");
1002 if (entity_dim == 0)
1005 auto& perm = _subentity_closure_perm_inv.at(entity_type);
1006 if (entity_dim == 1)
1008 if (entity_info & 1)
1013 else if (entity_dim == 2)
1016 if (entity_info & 1)
1022 for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
1029 throw std::runtime_error(
1030 "Invalid dimension for permute_subentity_closure");
1064 template <
typename T>
1065 void T_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1077 template <
typename T>
1078 void Tt_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1091 template <
typename T>
1092 void Tt_inv_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1104 template <
typename T>
1105 void Tinv_apply(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1117 template <
typename T>
1118 void Tt_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1129 template <
typename T>
1130 void T_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1142 template <
typename T>
1143 void Tinv_apply_right(std::span<T> u,
int n, std::uint32_t cell_info)
const;
1155 template <
typename T>
1164 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
points()
const
1220 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1231 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1234 return _dual_matrix;
1273 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
wcoeffs()
const
1283 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>, 4>&
1326 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>, 4>&
1335 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1353 return !_tensor_factors.empty();
1368 std::vector<std::vector<FiniteElement<F>>>
1372 throw std::runtime_error(
"Element has no tensor product representation.");
1373 return _tensor_factors;
1395 template <
typename T,
bool post>
1397 std::span<T> data,
int block_size, std::uint32_t cell_info,
1398 const std::map<
cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1401 using array2_t = std::pair<std::vector<F>, std::array<std::size_t, 2>>;
1402 using array3_t = std::pair<std::vector<F>, std::array<std::size_t, 3>>;
1404 = std::vector<std::pair<std::vector<std::size_t>, array2_t>>;
1412 template <
typename T,
bool post,
typename OP>
1414 transform_data(std::span<T> data,
int block_size, std::uint32_t cell_info,
1415 const std::map<cell::type, trans_data_t>& etrans, OP op)
const;
1424 std::size_t _cell_tdim;
1427 std::vector<std::vector<cell::type>> _cell_subentity_types;
1442 int _interpolation_nderivs;
1445 int _embedded_superdegree;
1448 int _embedded_subdegree;
1451 std::vector<std::size_t> _value_shape;
1464 std::pair<std::vector<F>, std::array<std::size_t, 2>> _coeffs;
1467 std::vector<std::vector<std::vector<int>>> _edofs;
1470 std::vector<std::vector<std::vector<int>>> _e_closure_dofs;
1473 std::map<cell::type, array3_t> _entity_transformations;
1480 std::pair<std::vector<F>, std::array<std::size_t, 2>> _points;
1484 std::array<std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>,
1489 std::pair<std::vector<F>, std::array<std::size_t, 2>> _matM;
1493 bool _dof_transformations_are_permutations;
1496 bool _dof_transformations_are_identity;
1501 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm;
1506 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm_inv;
1509 std::map<cell::type, trans_data_t> _etrans;
1512 std::map<cell::type, trans_data_t> _etransT;
1515 std::map<cell::type, trans_data_t> _etrans_inv;
1518 std::map<cell::type, trans_data_t> _etrans_invT;
1522 std::map<cell::type, std::vector<std::vector<std::size_t>>>
1523 _subentity_closure_perm;
1527 std::map<cell::type, std::vector<std::vector<std::size_t>>>
1528 _subentity_closure_perm_inv;
1532 bool _discontinuous;
1535 std::pair<std::vector<F>, std::array<std::size_t, 2>> _dual_matrix;
1543 std::vector<int> _dof_ordering;
1550 std::vector<std::vector<FiniteElement>> _tensor_factors;
1553 bool _interpolation_is_identity;
1557 std::pair<std::vector<F>, std::array<std::size_t, 2>> _wcoeffs;
1561 = std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>;
1562 std::array<array4_t, 4> _M;
1590 template <std::
floating_po
int T>
1592 cell::type cell_type,
const std::vector<std::size_t>& value_shape,
1593 impl::mdspan_t<const T, 2> wcoeffs,
1594 const std::array<std::vector<impl::mdspan_t<const T, 2>>, 4>& x,
1595 const std::array<std::vector<impl::mdspan_t<const T, 4>>, 4>& M,
1596 int interpolation_nderivs,
maps::type map_type,
1597 sobolev::space sobolev_space,
bool discontinuous,
int embedded_subdegree,
1611 template <std::
floating_po
int T>
1616 std::vector<int> dof_ordering = {});
1632 bool discontinuous);
1646 template <std::
floating_po
int T>
1647 std::vector<std::vector<FiniteElement<T>>>
1650 bool discontinuous, std::vector<int> dof_ordering);
1663 template <std::
floating_po
int T>
1674 template <std::
floating_po
int F>
1675 template <
typename T,
bool post>
1676 void FiniteElement<F>::permute_data(
1677 std::span<T> data,
int block_size, std::uint32_t cell_info,
1678 const std::map<
cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1681 if (_cell_tdim >= 2)
1685 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1689 auto& trans = eperm.at(cell::type::interval)[0];
1690 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1693 if (cell_info >> (face_start + e) & 1)
1701 if (_cell_tdim == 3)
1704 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1706 auto& trans = eperm.at(_cell_subentity_types[2][f]);
1709 if (!post and cell_info >> (3 * f) & 1)
1716 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1723 if (post and cell_info >> (3 * f) & 1)
1733 template <std::
floating_po
int F>
1734 template <
typename T,
bool post,
typename OP>
1735 void FiniteElement<F>::transform_data(
1736 std::span<T> data,
int block_size, std::uint32_t cell_info,
1737 const std::map<cell::type, trans_data_t>& etrans, OP op)
const
1739 if (_cell_tdim >= 2)
1743 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1745 for (
auto& edofs0 : _edofs[0])
1746 dofstart += edofs0.size();
1750 auto& [v_size_t, matrix] = etrans.at(cell::type::interval)[0];
1751 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1754 if (cell_info >> (face_start + e) & 1)
1756 op(std::span(v_size_t),
1757 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1758 dofstart, block_size);
1760 dofstart += _edofs[1][e].size();
1764 if (_cell_tdim == 3)
1767 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1769 auto& trans = etrans.at(_cell_subentity_types[2][f]);
1772 if (!post and cell_info >> (3 * f) & 1)
1774 const auto& m = trans[1];
1775 const auto& v_size_t = std::get<0>(m);
1776 const auto& matrix = std::get<1>(m);
1777 op(std::span(v_size_t),
1778 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1779 dofstart, block_size);
1783 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1785 const auto& m = trans[0];
1786 const auto& v_size_t = std::get<0>(m);
1787 const auto& matrix = std::get<1>(m);
1788 op(std::span(v_size_t),
1789 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1790 dofstart, block_size);
1794 if (post and cell_info >> (3 * f) & 1)
1796 const auto& m = trans[1];
1797 const auto& v_size_t = std::get<0>(m);
1798 const auto& matrix = std::get<1>(m);
1799 op(std::span(v_size_t),
1800 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1801 dofstart, block_size);
1804 dofstart += _edofs[2][f].size();
1810 template <std::
floating_po
int F>
1811 template <
typename T>
1813 std::uint32_t cell_info)
const
1815 if (_dof_transformations_are_identity)
1818 if (_dof_transformations_are_permutations)
1819 permute_data<T, false>(u, n, cell_info, _eperm);
1822 transform_data<T, false>(u, n, cell_info, _etrans,
1823 precompute::apply_matrix<F, T>);
1827 template <std::
floating_po
int F>
1828 template <
typename T>
1830 std::uint32_t cell_info)
const
1832 if (_dof_transformations_are_identity)
1834 else if (_dof_transformations_are_permutations)
1835 permute_data<T, true>(u, n, cell_info, _eperm_inv);
1838 transform_data<T, true>(u, n, cell_info, _etransT,
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, false>(u, n, cell_info, _eperm);
1854 transform_data<T, false>(u, n, cell_info, _etrans_invT,
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, true>(u, n, cell_info, _eperm_inv);
1870 transform_data<T, true>(u, n, cell_info, _etrans_inv,
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)
1884 assert(u.size() % n == 0);
1885 const int step = u.size() / n;
1886 for (
int i = 0; i < n; ++i)
1888 std::span<T> dblock(u.data() + i * step, step);
1889 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1894 transform_data<T, false>(u, n, cell_info, _etrans,
1895 precompute::apply_tranpose_matrix_right<F, T>);
1899 template <std::
floating_po
int F>
1900 template <
typename T>
1902 std::uint32_t cell_info)
const
1904 if (_dof_transformations_are_identity)
1906 else if (_dof_transformations_are_permutations)
1908 assert(u.size() % n == 0);
1909 const int step = u.size() / n;
1910 for (
int i = 0; i < n; ++i)
1912 std::span<T> dblock(u.data() + i * step, step);
1913 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1918 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1919 precompute::apply_tranpose_matrix_right<F, T>);
1923 template <std::
floating_po
int F>
1924 template <
typename T>
1926 std::uint32_t cell_info)
const
1928 if (_dof_transformations_are_identity)
1930 else if (_dof_transformations_are_permutations)
1932 assert(u.size() % n == 0);
1933 const int step = u.size() / n;
1934 for (
int i = 0; i < n; ++i)
1936 std::span<T> dblock(u.data() + i * step, step);
1937 permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1942 transform_data<T, true>(u, n, cell_info, _etransT,
1943 precompute::apply_tranpose_matrix_right<F, T>);
1947 template <std::
floating_po
int F>
1948 template <
typename T>
1950 std::uint32_t cell_info)
const
1952 if (_dof_transformations_are_identity)
1954 else if (_dof_transformations_are_permutations)
1956 assert(u.size() % n == 0);
1957 const int step = u.size() / n;
1958 for (
int i = 0; i < n; ++i)
1960 std::span<T> dblock(u.data() + i * step, step);
1961 permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1966 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1967 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:770
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:1901
const std::vector< int > & dof_ordering() const
Get dof layout.
Definition: finite-element.h:1386
std::vector< std::vector< FiniteElement< F > > > get_tensor_product_representation() const
Get the tensor product representation of this element.
Definition: finite-element.h:1369
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:1273
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:1327
std::pair< std::vector< F >, std::array< std::size_t, 3 > > base_transformations() const
Get the base transformations.
Definition: finite-element.cpp:1633
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:1812
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:1829
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:894
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & coefficient_matrix() const
Get the matrix of coefficients.
Definition: finite-element.h:1336
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & dual_matrix() const
Get the dual matrix.
Definition: finite-element.h:1232
bool dof_transformations_are_identity() const
Indicates is the dof transformations are all the identity.
Definition: finite-element.h:559
bool operator==(const FiniteElement &e) const
Check if two elements are the same.
Definition: finite-element.cpp:1452
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:990
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:507
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:1925
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:1351
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:676
int interpolation_nderivs() const
The number of derivatives needed when interpolating.
Definition: finite-element.h:1383
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:1541
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:502
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:625
int dim() const
Dimension of the finite element space.
Definition: finite-element.h:520
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:797
std::size_t hash() const
Get a unique hash of this element.
Definition: finite-element.cpp:1491
F scalar_type
Scalar type.
Definition: finite-element.h:145
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:828
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:1877
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:1845
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:857
element::lagrange_variant lagrange_variant() const
Lagrange variant of the element.
Definition: finite-element.h:528
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:1949
const std::vector< std::size_t > & value_shape() const
Element value tensor shape.
Definition: finite-element.h:514
int degree() const
Get the element polynomial degree.
Definition: finite-element.h:496
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:935
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & points() const
Return the interpolation points.
Definition: finite-element.h:1164
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:1861
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:365
FiniteElement & operator=(const FiniteElement &element)=default
Assignment operator.
sobolev::space sobolev_space() const
Underlying Sobolev space for this element.
Definition: finite-element.h:543
FiniteElement(const FiniteElement &element)=default
Copy constructor.
polyset::type polyset_type() const
Get the element polyset type.
Definition: finite-element.h:492
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:1221
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:1742
bool dof_transformations_are_permutations() const
Indicates if the degree-of-freedom transformations are all permutations.
Definition: finite-element.h:553
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:1702
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:662
cell::type cell_type() const
Get the element cell type.
Definition: finite-element.h:488
bool interpolation_is_identity() const
Indicates whether or not the interpolation matrix for this element is an identity matrix.
Definition: finite-element.h:1380
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:1284
FiniteElement & operator=(FiniteElement &&element)=default
Move assignment operator.
bool discontinuous() const
Indicates whether this element is the discontinuous variant.
Definition: finite-element.h:549
maps::type map_type() const
Map type for the element.
Definition: finite-element.h:539
element::dpc_variant dpc_variant() const
DPC variant of the element.
Definition: finite-element.h:535
element::family family() const
The finite element family.
Definition: finite-element.h:524
~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: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:522
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: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:195
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:401
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:351
std::string version()
Definition: finite-element.cpp:1780
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:332
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:611