Basix 0.10.0.0

Home     Installation     Demos     C++ docs     Python docs

finite-element.h
1 // Copyright (c) 2020-2024 Chris Richardson, Matthew Scroggs and Garth . Wells
2 // FEniCS Project
3 // SPDX-License-Identifier: MIT
4 
5 #pragma once
6 
7 #include "cell.h"
8 #include "element-families.h"
9 #include "maps.h"
10 #include "mdspan.hpp"
11 #include "polyset.h"
12 #include "precompute.h"
13 #include "sobolev-spaces.h"
14 #include <array>
15 #include <concepts>
16 #include <cstdint>
17 #include <functional>
18 #include <map>
19 #include <numeric>
20 #include <span>
21 #include <string>
22 #include <tuple>
23 #include <utility>
24 #include <vector>
25 
27 namespace basix
28 {
29 namespace impl
30 {
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>
35 using mdarray_t
36  = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::mdarray<
37  T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
38 
41 template <typename T>
42 std::array<std::vector<mdspan_t<const T, 2>>, 4>
43 to_mdspan(std::array<std::vector<mdarray_t<T, 2>>, 4>& x)
44 {
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());
49 
50  return x1;
51 }
52 
55 template <typename T>
56 std::array<std::vector<mdspan_t<const T, 4>>, 4>
57 to_mdspan(std::array<std::vector<mdarray_t<T, 4>>, 4>& M)
58 {
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());
63 
64  return M1;
65 }
66 
69 template <typename T>
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)
73 {
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]));
78 
79  return x1;
80 }
81 
84 template <typename T>
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)
88 {
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]));
93 
94  return M1;
95 }
96 
97 } // namespace impl
98 
99 namespace element
100 {
102 template <typename T, std::size_t d>
103 using mdspan_t = impl::mdspan_t<T, d>;
104 
120 template <std::floating_point 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>>
125 make_discontinuous(const std::array<std::vector<mdspan_t<const T, 2>>, 4>& x,
126  const std::array<std::vector<mdspan_t<const T, 4>>, 4>& M,
127  std::size_t tdim, std::size_t value_size);
128 
129 } // namespace element
130 
136 template <std::floating_point F>
138 {
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>>;
142 
143 public:
145  using scalar_type = F;
146 
320  polyset::type poly_type, int degree,
321  const std::vector<std::size_t>& value_shape,
322  mdspan_t<const F, 2> wcoeffs,
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,
328  element::lagrange_variant lvariant,
329  element::dpc_variant dvariant,
330  std::vector<int> dof_ordering = {});
331 
333  FiniteElement(const FiniteElement& element) = default;
334 
336  FiniteElement(FiniteElement&& element) = default;
337 
339  ~FiniteElement() = default;
340 
342  FiniteElement& operator=(const FiniteElement& element) = default;
343 
345  FiniteElement& operator=(FiniteElement&& element) = default;
346 
351  bool operator==(const FiniteElement& e) const;
352 
354  std::size_t hash() const;
355 
365  std::array<std::size_t, 4> tabulate_shape(std::size_t nd,
366  std::size_t num_points) const
367  {
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)
372  ndsize /= 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};
377  }
378 
400  std::pair<std::vector<F>, std::array<std::size_t, 4>>
401  tabulate(int nd, impl::mdspan_t<const F, 2> x) const;
402 
426  std::pair<std::vector<F>, std::array<std::size_t, 4>>
427  tabulate(int nd, std::span<const F> x,
428  std::array<std::size_t, 2> shape) const;
429 
455  void tabulate(int nd, impl::mdspan_t<const F, 2> x,
456  mdspan_t<F, 4> basis) const;
457 
483  void tabulate(int nd, std::span<const F> x, std::array<std::size_t, 2> xshape,
484  std::span<F> basis) const;
485 
488  cell::type cell_type() const { return _cell_type; }
489 
492  polyset::type polyset_type() const { return _poly_type; }
493 
496  int degree() const { return _degree; }
497 
502  int embedded_superdegree() const { return _embedded_superdegree; }
503 
507  int embedded_subdegree() const { return _embedded_subdegree; }
508 
514  const std::vector<std::size_t>& value_shape() const { return _value_shape; }
515 
520  int dim() const { return _coeffs.second[0]; }
521 
524  element::family family() const { return _family; }
525 
529  {
530  return _lagrange_variant;
531  }
532 
535  element::dpc_variant dpc_variant() const { return _dpc_variant; }
536 
539  maps::type map_type() const { return _map_type; }
540 
543  sobolev::space sobolev_space() const { return _sobolev_space; }
544 
549  bool discontinuous() const { return _discontinuous; }
550 
554  {
555  return _dof_transformations_are_permutations;
556  }
557 
560  {
561  return _dof_transformations_are_identity;
562  }
563 
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;
581 
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;
592 
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
626  {
627  switch (_map_type)
628  {
629  case maps::type::identity:
630  return [](O& u, const P& U, const Q&, F, const R&)
631  {
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)
636  u(i, j) = U(i, j);
637  };
638  case maps::type::covariantPiola:
639  return [](O& u, const P& U, const Q& J, F detJ, const R& K)
640  { maps::covariant_piola(u, U, J, detJ, K); };
641  case maps::type::contravariantPiola:
642  return [](O& u, const P& U, const Q& J, F detJ, const R& K)
643  { maps::contravariant_piola(u, U, J, detJ, K); };
644  case maps::type::doubleCovariantPiola:
645  return [](O& u, const P& U, const Q& J, F detJ, const R& K)
646  { maps::double_covariant_piola(u, U, J, detJ, K); };
647  case maps::type::doubleContravariantPiola:
648  return [](O& u, const P& U, const Q& J, F detJ, const R& K)
649  { maps::double_contravariant_piola(u, U, J, detJ, K); };
650  default:
651  throw std::runtime_error("Map not implemented");
652  }
653  }
654 
662  const std::vector<std::vector<std::vector<int>>>& entity_dofs() const
663  {
664  return _edofs;
665  }
666 
676  const std::vector<std::vector<std::vector<int>>>& entity_closure_dofs() const
677  {
678  return _e_closure_dofs;
679  }
680 
762  std::pair<std::vector<F>, std::array<std::size_t, 3>>
763  base_transformations() const;
764 
769  std::map<cell::type, std::pair<std::vector<F>, std::array<std::size_t, 3>>>
771  {
772  return _entity_transformations;
773  }
774 
797  void permute(std::span<std::int32_t> d, std::uint32_t cell_info) const
798  {
799  if (!_dof_transformations_are_permutations)
800  {
801  throw std::runtime_error(
802  "The DOF transformations for this element are not permutations");
803  }
804 
805  if (_dof_transformations_are_identity)
806  return;
807  else
808  permute_data<std::int32_t, false>(d, 1, cell_info, _eperm);
809  }
810 
828  void permute_inv(std::span<std::int32_t> d, std::uint32_t cell_info) const
829  {
830  if (!_dof_transformations_are_permutations)
831  {
832  throw std::runtime_error(
833  "The DOF transformations for this element are not permutations");
834  }
835 
836  if (_dof_transformations_are_identity)
837  return;
838  else
839  permute_data<std::int32_t, true>(d, 1, cell_info, _eperm_inv);
840  }
841 
857  void permute_subentity_closure(std::span<std::int32_t> d,
858  std::uint32_t cell_info,
859  cell::type entity_type, int entity_index) const
860  {
861  const int entity_dim = cell::topological_dimension(entity_type);
862 
863  int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
864 
865  std::uint32_t entity_info = 0;
866  switch (entity_dim)
867  {
868  case 0:
869  return;
870  case 1:
871  entity_info = cell_info >> (face_start + entity_index) & 1;
872  break;
873  case 2:
874  entity_info = cell_info >> (3 * entity_index) & 7;
875  break;
876  default:
877  throw std::runtime_error("Unsupported cell dimension");
878  }
879  permute_subentity_closure(d, entity_info, entity_type);
880  }
881 
894  void permute_subentity_closure_inv(std::span<std::int32_t> d,
895  std::uint32_t cell_info,
896  cell::type entity_type,
897  int entity_index) const
898  {
899  const int entity_dim = cell::topological_dimension(entity_type);
900 
901  int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
902 
903  std::uint32_t entity_info;
904  switch (entity_dim)
905  {
906  case 0:
907  return;
908  case 1:
909  entity_info = cell_info >> (face_start + entity_index) & 1;
910  break;
911  case 2:
912  entity_info = cell_info >> (3 * entity_index) & 7;
913  break;
914  default:
915  throw std::runtime_error("Unsupported cell dimension");
916  }
917  permute_subentity_closure_inv(d, entity_info, entity_type);
918  }
919 
934 
935  void permute_subentity_closure(std::span<std::int32_t> d,
936  std::uint32_t entity_info,
937  cell::type entity_type) const
938  {
939  if (!_dof_transformations_are_permutations)
940  {
941  throw std::runtime_error(
942  "The DOF transformations for this element are not permutations");
943  }
944 
945  const int entity_dim = cell::topological_dimension(entity_type);
946 
947  if (entity_dim == 0)
948  return;
949 
950  auto& perm = _subentity_closure_perm.at(entity_type);
951  if (entity_dim == 1)
952  {
953  if (entity_info & 1)
954  {
955  precompute::apply_permutation(perm[0], d);
956  }
957  }
958  else if (entity_dim == 2)
959  {
960  // Rotate a face
961  for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
962  {
963  precompute::apply_permutation(perm[0], d);
964  }
965 
966  // Reflect a face (post rotate)
967  if (entity_info & 1)
968  {
969  precompute::apply_permutation(perm[1], d);
970  }
971  }
972  else
973  {
974  throw std::runtime_error(
975  "Invalid dimension for permute_subentity_closure");
976  }
977  }
978 
990  void permute_subentity_closure_inv(std::span<std::int32_t> d,
991  std::uint32_t entity_info,
992  cell::type entity_type) const
993  {
994  if (!_dof_transformations_are_permutations)
995  {
996  throw std::runtime_error(
997  "The DOF transformations for this element are not permutations");
998  }
999 
1000  const int entity_dim = cell::topological_dimension(entity_type);
1001 
1002  if (entity_dim == 0)
1003  return;
1004 
1005  auto& perm = _subentity_closure_perm_inv.at(entity_type);
1006  if (entity_dim == 1)
1007  {
1008  if (entity_info & 1)
1009  {
1010  precompute::apply_permutation(perm[0], d);
1011  }
1012  }
1013  else if (entity_dim == 2)
1014  {
1015  // Reflect a face (pre rotate)
1016  if (entity_info & 1)
1017  {
1018  precompute::apply_permutation(perm[1], d);
1019  }
1020 
1021  // Rotate a face
1022  for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
1023  {
1024  precompute::apply_permutation(perm[0], d);
1025  }
1026  }
1027  else
1028  {
1029  throw std::runtime_error(
1030  "Invalid dimension for permute_subentity_closure");
1031  }
1032  }
1033 
1064  template <typename T>
1065  void T_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
1066 
1077  template <typename T>
1078  void Tt_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
1079 
1091  template <typename T>
1092  void Tt_inv_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
1093 
1104  template <typename T>
1105  void Tinv_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
1106 
1117  template <typename T>
1118  void Tt_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
1119 
1129  template <typename T>
1130  void T_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
1131 
1142  template <typename T>
1143  void Tinv_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
1144 
1155  template <typename T>
1156  void Tt_inv_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
1157 
1164  const std::pair<std::vector<F>, std::array<std::size_t, 2>>& points() const
1165  {
1166  return _points;
1167  }
1168 
1220  const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1222  {
1223  return _matM;
1224  }
1225 
1231  const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1232  dual_matrix() const
1233  {
1234  return _dual_matrix;
1235  }
1236 
1273  const std::pair<std::vector<F>, std::array<std::size_t, 2>>& wcoeffs() const
1274  {
1275  return _wcoeffs;
1276  }
1277 
1282  const std::array<
1283  std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>, 4>&
1284  x() const
1285  {
1286  return _x;
1287  }
1288 
1325  const std::array<
1326  std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>, 4>&
1327  M() const
1328  {
1329  return _M;
1330  }
1331 
1335  const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1337  {
1338  return _coeffs;
1339  }
1340 
1352  {
1353  return !_tensor_factors.empty();
1354  }
1355 
1368  std::vector<std::vector<FiniteElement<F>>>
1370  {
1372  throw std::runtime_error("Element has no tensor product representation.");
1373  return _tensor_factors;
1374  }
1375 
1380  bool interpolation_is_identity() const { return _interpolation_is_identity; }
1381 
1383  int interpolation_nderivs() const { return _interpolation_nderivs; }
1384 
1386  const std::vector<int>& dof_ordering() const { return _dof_ordering; }
1387 
1388 private:
1395  template <typename T, bool post>
1396  void permute_data(
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)
1399  const;
1400 
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>>;
1403  using trans_data_t
1404  = std::vector<std::pair<std::vector<std::size_t>, array2_t>>;
1405 
1412  template <typename T, bool post, typename OP>
1413  void
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;
1416 
1417  // Cell type
1418  cell::type _cell_type;
1419 
1420  // Polyset type
1421  polyset::type _poly_type;
1422 
1423  // Topological dimension of the cell
1424  std::size_t _cell_tdim;
1425 
1426  // Topological dimension of the cell
1427  std::vector<std::vector<cell::type>> _cell_subentity_types;
1428 
1429  // Finite element family
1430  element::family _family;
1431 
1432  // Lagrange variant
1433  element::lagrange_variant _lagrange_variant;
1434 
1435  // DPC variant
1436  element::dpc_variant _dpc_variant;
1437 
1438  // Degree that was input when creating the element
1439  int _degree;
1440 
1441  // Degree
1442  int _interpolation_nderivs;
1443 
1444  // Highest degree polynomial in element's polyset
1445  int _embedded_superdegree;
1446 
1447  // Highest degree space that is a subspace of element's polyset
1448  int _embedded_subdegree;
1449 
1450  // Value shape
1451  std::vector<std::size_t> _value_shape;
1452 
1454  maps::type _map_type;
1455 
1457  sobolev::space _sobolev_space;
1458 
1459  // Shape function coefficient of expansion sets on cell. If shape
1460  // function is given by @f$\psi_i = \sum_{k} \phi_{k}
1461  // \alpha^{i}_{k}@f$, then _coeffs(i, j) = @f$\alpha^i_k@f$. ie
1462  // _coeffs.row(i) are the expansion coefficients for shape function i
1463  // (@f$\psi_{i}@f$).
1464  std::pair<std::vector<F>, std::array<std::size_t, 2>> _coeffs;
1465 
1466  // Dofs associated with each cell (sub-)entity
1467  std::vector<std::vector<std::vector<int>>> _edofs;
1468 
1469  // Dofs associated with the closdure of each cell (sub-)entity
1470  std::vector<std::vector<std::vector<int>>> _e_closure_dofs;
1471 
1472  // Entity transformations
1473  std::map<cell::type, array3_t> _entity_transformations;
1474 
1475  // Set of points used for point evaluation
1476  // Experimental - currently used for an implementation of
1477  // "tabulate_dof_coordinates" Most useful for Lagrange. This may change or go
1478  // away. For non-Lagrange elements, these points will be used in combination
1479  // with _interpolation_matrix to perform interpolation
1480  std::pair<std::vector<F>, std::array<std::size_t, 2>> _points;
1481 
1482  // Interpolation points on the cell. The shape is (entity_dim, num
1483  // entities of given dimension, num_points, tdim)
1484  std::array<std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>,
1485  4>
1486  _x;
1487 
1489  std::pair<std::vector<F>, std::array<std::size_t, 2>> _matM;
1490 
1491  // Indicates whether or not the DOF transformations are all
1492  // permutations
1493  bool _dof_transformations_are_permutations;
1494 
1495  // Indicates whether or not the DOF transformations are all identity
1496  bool _dof_transformations_are_identity;
1497 
1498  // The entity permutations (factorised). This will only be set if
1499  // _dof_transformations_are_permutations is True and
1500  // _dof_transformations_are_identity is False
1501  std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm;
1502 
1503  // The reverse entity permutations (factorised). This will only be set
1504  // if _dof_transformations_are_permutations is True and
1505  // _dof_transformations_are_identity is False
1506  std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm_inv;
1507 
1508  // The entity transformations in precomputed form
1509  std::map<cell::type, trans_data_t> _etrans;
1510 
1511  // The transposed entity transformations in precomputed form
1512  std::map<cell::type, trans_data_t> _etransT;
1513 
1514  // The inverse entity transformations in precomputed form
1515  std::map<cell::type, trans_data_t> _etrans_inv;
1516 
1517  // The inverse transpose entity transformations in precomputed form
1518  std::map<cell::type, trans_data_t> _etrans_invT;
1519 
1520  // The subentity closure permutations (factorised). This will only be set if
1521  // _dof_transformations_are_permutations is True
1522  std::map<cell::type, std::vector<std::vector<std::size_t>>>
1523  _subentity_closure_perm;
1524 
1525  // The inverse subentity closure permutations (factorised). This will only be
1526  // set if _dof_transformations_are_permutations is True
1527  std::map<cell::type, std::vector<std::vector<std::size_t>>>
1528  _subentity_closure_perm_inv;
1529 
1530  // Indicates whether or not this is the discontinuous version of the
1531  // element
1532  bool _discontinuous;
1533 
1534  // The dual matrix
1535  std::pair<std::vector<F>, std::array<std::size_t, 2>> _dual_matrix;
1536 
1537  // Dof reordering for different element dof layout compatibility.
1538  // The reference basix layout is ordered by entity, i.e. dofs on
1539  // vertices, followed by edges, faces, then internal dofs.
1540  // _dof_ordering stores the map to the new order required, e.g.
1541  // for a P2 triangle, _dof_ordering=[0 3 5 1 2 4] will place
1542  // dofs 0, 3, 5 on the vertices and 1, 2, 4, on the edges.
1543  std::vector<int> _dof_ordering;
1544 
1545  // Tensor product representation
1546  // Entries of tuple are (list of elements on an interval, permutation
1547  // of DOF numbers)
1548  // @todo: For vector-valued elements, a tensor product type and a
1549  // scaling factor may additionally be needed.
1550  std::vector<std::vector<FiniteElement>> _tensor_factors;
1551 
1552  // Is the interpolation matrix an identity?
1553  bool _interpolation_is_identity;
1554 
1555  // The coefficients that define the polynomial set in terms of the
1556  // orthonormal polynomials
1557  std::pair<std::vector<F>, std::array<std::size_t, 2>> _wcoeffs;
1558 
1559  // Interpolation matrices for each entity
1560  using array4_t
1561  = std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>;
1562  std::array<array4_t, 4> _M;
1563 };
1564 
1590 template <std::floating_point T>
1591 FiniteElement<T> create_custom_element(
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,
1598  int embedded_superdegree, polyset::type poly_type);
1599 
1611 template <std::floating_point T>
1612 FiniteElement<T> create_element(element::family family, cell::type cell,
1613  int degree, element::lagrange_variant lvariant,
1614  element::dpc_variant dvariant,
1615  bool discontinuous,
1616  std::vector<int> dof_ordering = {});
1617 
1629 std::vector<int> tp_dof_ordering(element::family family, cell::type cell,
1630  int degree, element::lagrange_variant lvariant,
1631  element::dpc_variant dvariant,
1632  bool discontinuous);
1633 
1646 template <std::floating_point T>
1647 std::vector<std::vector<FiniteElement<T>>>
1648 tp_factors(element::family family, cell::type cell, int degree,
1650  bool discontinuous, std::vector<int> dof_ordering);
1651 
1663 template <std::floating_point T>
1664 FiniteElement<T>
1665 create_tp_element(element::family family, cell::type cell, int degree,
1666  element::lagrange_variant lvariant,
1667  element::dpc_variant dvariant, bool discontinuous);
1668 
1671 std::string version();
1672 
1673 //-----------------------------------------------------------------------------
1674 template <std::floating_point 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)
1679  const
1680 {
1681  if (_cell_tdim >= 2)
1682  {
1683  // This assumes 3 bits are used per face. This will need updating if 3D
1684  // cells with faces with more than 4 sides are implemented
1685  int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1686 
1687  // Permute DOFs on edges
1688  {
1689  auto& trans = eperm.at(cell::type::interval)[0];
1690  for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1691  {
1692  // Reverse an edge
1693  if (cell_info >> (face_start + e) & 1)
1694  {
1695  precompute::apply_permutation_mapped(trans, data, _edofs[1][e],
1696  block_size);
1697  }
1698  }
1699  }
1700 
1701  if (_cell_tdim == 3)
1702  {
1703  // Permute DOFs on faces
1704  for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1705  {
1706  auto& trans = eperm.at(_cell_subentity_types[2][f]);
1707 
1708  // Reflect a face (pre rotate)
1709  if (!post and cell_info >> (3 * f) & 1)
1710  {
1711  precompute::apply_permutation_mapped(trans[1], data, _edofs[2][f],
1712  block_size);
1713  }
1714 
1715  // Rotate a face
1716  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1717  {
1718  precompute::apply_permutation_mapped(trans[0], data, _edofs[2][f],
1719  block_size);
1720  }
1721 
1722  // Reflect a face (post rotate)
1723  if (post and cell_info >> (3 * f) & 1)
1724  {
1725  precompute::apply_permutation_mapped(trans[1], data, _edofs[2][f],
1726  block_size);
1727  }
1728  }
1729  }
1730  }
1731 }
1732 //-----------------------------------------------------------------------------
1733 template <std::floating_point 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
1738 {
1739  if (_cell_tdim >= 2)
1740  {
1741  // This assumes 3 bits are used per face. This will need updating if
1742  // 3D cells with faces with more than 4 sides are implemented
1743  int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1744  int dofstart = 0;
1745  for (auto& edofs0 : _edofs[0])
1746  dofstart += edofs0.size();
1747 
1748  // Transform DOFs on edges
1749  {
1750  auto& [v_size_t, matrix] = etrans.at(cell::type::interval)[0];
1751  for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1752  {
1753  // Reverse an edge
1754  if (cell_info >> (face_start + e) & 1)
1755  {
1756  op(std::span(v_size_t),
1757  mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1758  dofstart, block_size);
1759  }
1760  dofstart += _edofs[1][e].size();
1761  }
1762  }
1763 
1764  if (_cell_tdim == 3)
1765  {
1766  // Permute DOFs on faces
1767  for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1768  {
1769  auto& trans = etrans.at(_cell_subentity_types[2][f]);
1770 
1771  // Reflect a face (pre rotation)
1772  if (!post and cell_info >> (3 * f) & 1)
1773  {
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);
1780  }
1781 
1782  // Rotate a face
1783  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1784  {
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);
1791  }
1792 
1793  // Reflect a face (post rotation)
1794  if (post and cell_info >> (3 * f) & 1)
1795  {
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);
1802  }
1803 
1804  dofstart += _edofs[2][f].size();
1805  }
1806  }
1807  }
1808 }
1809 //-----------------------------------------------------------------------------
1810 template <std::floating_point F>
1811 template <typename T>
1812 void FiniteElement<F>::T_apply(std::span<T> u, int n,
1813  std::uint32_t cell_info) const
1814 {
1815  if (_dof_transformations_are_identity)
1816  return;
1817 
1818  if (_dof_transformations_are_permutations)
1819  permute_data<T, false>(u, n, cell_info, _eperm);
1820  else
1821  {
1822  transform_data<T, false>(u, n, cell_info, _etrans,
1823  precompute::apply_matrix<F, T>);
1824  }
1825 }
1826 //-----------------------------------------------------------------------------
1827 template <std::floating_point F>
1828 template <typename T>
1829 void FiniteElement<F>::Tt_apply(std::span<T> u, int n,
1830  std::uint32_t cell_info) const
1831 {
1832  if (_dof_transformations_are_identity)
1833  return;
1834  else if (_dof_transformations_are_permutations)
1835  permute_data<T, true>(u, n, cell_info, _eperm_inv);
1836  else
1837  {
1838  transform_data<T, true>(u, n, cell_info, _etransT,
1839  precompute::apply_matrix<F, T>);
1840  }
1841 }
1842 //-----------------------------------------------------------------------------
1843 template <std::floating_point F>
1844 template <typename T>
1845 void FiniteElement<F>::Tt_inv_apply(std::span<T> u, int n,
1846  std::uint32_t cell_info) const
1847 {
1848  if (_dof_transformations_are_identity)
1849  return;
1850  else if (_dof_transformations_are_permutations)
1851  permute_data<T, false>(u, n, cell_info, _eperm);
1852  else
1853  {
1854  transform_data<T, false>(u, n, cell_info, _etrans_invT,
1855  precompute::apply_matrix<F, T>);
1856  }
1857 }
1858 //-----------------------------------------------------------------------------
1859 template <std::floating_point F>
1860 template <typename T>
1861 void FiniteElement<F>::Tinv_apply(std::span<T> u, int n,
1862  std::uint32_t cell_info) const
1863 {
1864  if (_dof_transformations_are_identity)
1865  return;
1866  else if (_dof_transformations_are_permutations)
1867  permute_data<T, true>(u, n, cell_info, _eperm_inv);
1868  else
1869  {
1870  transform_data<T, true>(u, n, cell_info, _etrans_inv,
1871  precompute::apply_matrix<F, T>);
1872  }
1873 }
1874 //-----------------------------------------------------------------------------
1875 template <std::floating_point F>
1876 template <typename T>
1877 void FiniteElement<F>::Tt_apply_right(std::span<T> u, int n,
1878  std::uint32_t cell_info) const
1879 {
1880  if (_dof_transformations_are_identity)
1881  return;
1882  else if (_dof_transformations_are_permutations)
1883  {
1884  assert(u.size() % n == 0);
1885  const int step = u.size() / n;
1886  for (int i = 0; i < n; ++i)
1887  {
1888  std::span<T> dblock(u.data() + i * step, step);
1889  permute_data<T, false>(dblock, 1, cell_info, _eperm);
1890  }
1891  }
1892  else
1893  {
1894  transform_data<T, false>(u, n, cell_info, _etrans,
1895  precompute::apply_tranpose_matrix_right<F, T>);
1896  }
1897 }
1898 //-----------------------------------------------------------------------------
1899 template <std::floating_point F>
1900 template <typename T>
1901 void FiniteElement<F>::Tinv_apply_right(std::span<T> u, int n,
1902  std::uint32_t cell_info) const
1903 {
1904  if (_dof_transformations_are_identity)
1905  return;
1906  else if (_dof_transformations_are_permutations)
1907  {
1908  assert(u.size() % n == 0);
1909  const int step = u.size() / n;
1910  for (int i = 0; i < n; ++i)
1911  {
1912  std::span<T> dblock(u.data() + i * step, step);
1913  permute_data<T, false>(dblock, 1, cell_info, _eperm);
1914  }
1915  }
1916  else
1917  {
1918  transform_data<T, false>(u, n, cell_info, _etrans_invT,
1919  precompute::apply_tranpose_matrix_right<F, T>);
1920  }
1921 }
1922 //-----------------------------------------------------------------------------
1923 template <std::floating_point F>
1924 template <typename T>
1925 void FiniteElement<F>::T_apply_right(std::span<T> u, int n,
1926  std::uint32_t cell_info) const
1927 {
1928  if (_dof_transformations_are_identity)
1929  return;
1930  else if (_dof_transformations_are_permutations)
1931  {
1932  assert(u.size() % n == 0);
1933  const int step = u.size() / n;
1934  for (int i = 0; i < n; ++i)
1935  {
1936  std::span<T> dblock(u.data() + i * step, step);
1937  permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1938  }
1939  }
1940  else
1941  {
1942  transform_data<T, true>(u, n, cell_info, _etransT,
1943  precompute::apply_tranpose_matrix_right<F, T>);
1944  }
1945 }
1946 //-----------------------------------------------------------------------------
1947 template <std::floating_point F>
1948 template <typename T>
1949 void FiniteElement<F>::Tt_inv_apply_right(std::span<T> u, int n,
1950  std::uint32_t cell_info) const
1951 {
1952  if (_dof_transformations_are_identity)
1953  return;
1954  else if (_dof_transformations_are_permutations)
1955  {
1956  assert(u.size() % n == 0);
1957  const int step = u.size() / n;
1958  for (int i = 0; i < n; ++i)
1959  {
1960  std::span<T> dblock(u.data() + i * step, step);
1961  permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1962  }
1963  }
1964  else
1965  {
1966  transform_data<T, true>(u, n, cell_info, _etrans_inv,
1967  precompute::apply_tranpose_matrix_right<F, T>);
1968  }
1969 }
1970 //-----------------------------------------------------------------------------
1971 
1972 } // namespace basix
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