Basix 0.4.1

Home     Installation     Demos     C++ docs     Python docs

finite-element.h
1 // Copyright (c) 2020 Chris Richardson
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 "precompute.h"
11 #include <array>
12 #include <functional>
13 #include <numeric>
14 #include <string>
15 #include <vector>
16 #include <xtensor/xtensor.hpp>
17 #include <xtensor/xview.hpp>
18 #include <xtl/xspan.hpp>
19 
21 namespace basix
22 {
23 
24 namespace element
25 {
26 
39 std::tuple<std::array<std::vector<xt::xtensor<double, 2>>, 4>,
40  std::array<std::vector<xt::xtensor<double, 3>>, 4>>
41 make_discontinuous(const std::array<std::vector<xt::xtensor<double, 2>>, 4>& x,
42  const std::array<std::vector<xt::xtensor<double, 3>>, 4>& M,
43  int tdim, int value_size);
44 
45 } // namespace element
46 
48 
53 {
54 
55 public:
223  const std::vector<std::size_t>& value_shape,
224  const xt::xtensor<double, 2>& wcoeffs,
225  const std::array<std::vector<xt::xtensor<double, 2>>, 4>& x,
226  const std::array<std::vector<xt::xtensor<double, 3>>, 4>& M,
229  element::dpc_variant dvariant,
230  std::vector<std::tuple<std::vector<FiniteElement>, std::vector<int>>>
231  tensor_factors
232  = {});
233 
237  const std::vector<std::size_t>& value_shape,
238  const xt::xtensor<double, 2>& wcoeffs,
239  const std::array<std::vector<xt::xtensor<double, 2>>, 4>& x,
240  const std::array<std::vector<xt::xtensor<double, 3>>, 4>& M,
243  std::vector<std::tuple<std::vector<FiniteElement>, std::vector<int>>>
244  tensor_factors
245  = {});
246 
250  const std::vector<std::size_t>& value_shape,
251  const xt::xtensor<double, 2>& wcoeffs,
252  const std::array<std::vector<xt::xtensor<double, 2>>, 4>& x,
253  const std::array<std::vector<xt::xtensor<double, 3>>, 4>& M,
255  int highest_degree, element::dpc_variant dvariant,
256  std::vector<std::tuple<std::vector<FiniteElement>, std::vector<int>>>
257  tensor_factors
258  = {});
259 
263  const std::vector<std::size_t>& value_shape,
264  const xt::xtensor<double, 2>& wcoeffs,
265  const std::array<std::vector<xt::xtensor<double, 2>>, 4>& x,
266  const std::array<std::vector<xt::xtensor<double, 3>>, 4>& M,
268  int highest_degree,
269  std::vector<std::tuple<std::vector<FiniteElement>, std::vector<int>>>
270  tensor_factors
271  = {});
272 
274  FiniteElement(const FiniteElement& element) = default;
275 
277  FiniteElement(FiniteElement&& element) = default;
278 
280  ~FiniteElement() = default;
281 
283  FiniteElement& operator=(const FiniteElement& element) = default;
284 
286  FiniteElement& operator=(FiniteElement&& element) = default;
287 
292  bool operator==(const FiniteElement& e) const;
293 
303  std::array<std::size_t, 4> tabulate_shape(std::size_t nd,
304  std::size_t num_points) const;
305 
327  xt::xtensor<double, 4> tabulate(int nd,
328  const xt::xtensor<double, 2>& x) const;
329 
355  void tabulate(int nd, const xt::xtensor<double, 2>& x,
356  xt::xtensor<double, 4>& basis) const;
357 
360  cell::type cell_type() const;
361 
364  int degree() const;
365 
370  int highest_degree() const;
371 
375  int highest_complete_degree() const;
376 
380  const std::vector<int>& value_shape() const;
381 
385  int dim() const;
386 
389  element::family family() const;
390 
394 
398 
401  maps::type map_type() const;
402 
406  bool discontinuous() const;
407 
411 
415 
429  xt::xtensor<double, 3> push_forward(const xt::xtensor<double, 3>& U,
430  const xt::xtensor<double, 3>& J,
431  const xtl::span<const double>& detJ,
432  const xt::xtensor<double, 3>& K) const;
433 
440  xt::xtensor<double, 3> pull_back(const xt::xtensor<double, 3>& u,
441  const xt::xtensor<double, 3>& J,
442  const xtl::span<const double>& detJ,
443  const xt::xtensor<double, 3>& K) const;
444 
476  template <typename O, typename P, typename Q, typename R>
477  std::function<void(O&, const P&, const Q&, double, const R&)> map_fn() const
478  {
479  switch (_map_type)
480  {
481  case maps::type::identity:
482  return [](O& u, const P& U, const Q&, double, const R&) { u.assign(U); };
483  case maps::type::L2Piola:
484  return [](O& u, const P& U, const Q& J, double detJ, const R& K)
485  { maps::l2_piola(u, U, J, detJ, K); };
486  case maps::type::covariantPiola:
487  return [](O& u, const P& U, const Q& J, double detJ, const R& K)
488  { maps::covariant_piola(u, U, J, detJ, K); };
489  case maps::type::contravariantPiola:
490  return [](O& u, const P& U, const Q& J, double detJ, const R& K)
491  { maps::contravariant_piola(u, U, J, detJ, K); };
492  case maps::type::doubleCovariantPiola:
493  return [](O& u, const P& U, const Q& J, double detJ, const R& K)
494  { maps::double_covariant_piola(u, U, J, detJ, K); };
495  case maps::type::doubleContravariantPiola:
496  return [](O& u, const P& U, const Q& J, double detJ, const R& K)
497  { maps::double_contravariant_piola(u, U, J, detJ, K); };
498  default:
499  throw std::runtime_error("Map not implemented");
500  }
501  }
502 
515  const std::vector<std::vector<int>>& num_entity_dofs() const;
516 
522  const std::vector<std::vector<int>>& num_entity_closure_dofs() const;
523 
530  const std::vector<std::vector<std::vector<int>>>& entity_dofs() const;
531 
538  const std::vector<std::vector<std::vector<int>>>& entity_closure_dofs() const;
539 
619  xt::xtensor<double, 3> base_transformations() const;
620 
623  std::map<cell::type, xt::xtensor<double, 3>> entity_transformations() const;
624 
632  void permute_dofs(const xtl::span<std::int32_t>& dofs,
633  std::uint32_t cell_info) const;
634 
642  void unpermute_dofs(const xtl::span<std::int32_t>& dofs,
643  std::uint32_t cell_info) const;
644 
653  template <typename T>
654  void apply_dof_transformation(const xtl::span<T>& data, int block_size,
655  std::uint32_t cell_info) const;
656 
665  template <typename T>
666  void apply_transpose_dof_transformation(const xtl::span<T>& data,
667  int block_size,
668  std::uint32_t cell_info) const;
669 
678  template <typename T>
680  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const;
681 
690  template <typename T>
691  void apply_inverse_dof_transformation(const xtl::span<T>& data,
692  int block_size,
693  std::uint32_t cell_info) const;
694 
703  template <typename T>
704  void apply_dof_transformation_to_transpose(const xtl::span<T>& data,
705  int block_size,
706  std::uint32_t cell_info) const;
707 
716  template <typename T>
718  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const;
719 
728  template <typename T>
730  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const;
731 
740  template <typename T>
742  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const;
743 
748  const xt::xtensor<double, 2>& points() const;
749 
799  const xt::xtensor<double, 2>& interpolation_matrix() const;
800 
806  const xt::xtensor<double, 2>& dual_matrix() const;
807 
839  const xt::xtensor<double, 2>& wcoeffs() const;
840 
844  const std::array<std::vector<xt::xtensor<double, 2>>, 4>& x() const;
845 
878  const std::array<std::vector<xt::xtensor<double, 3>>, 4>& M() const;
879 
885  const xt::xtensor<double, 2>& coefficient_matrix() const;
886 
887 
891 
903  std::vector<std::tuple<std::vector<FiniteElement>, std::vector<int>>>
905 
908  bool interpolation_is_identity() const;
909 
910 private:
911  // Cell type
912  cell::type _cell_type;
913 
914  // Topological dimension of the cell
915  std::size_t _cell_tdim;
916 
917  // Topological dimension of the cell
918  std::vector<std::vector<cell::type>> _cell_subentity_types;
919 
920  // Finite element family
921  element::family _family;
922 
923  // Lagrange variant
924  element::lagrange_variant _lagrange_variant;
925 
926  // Lagrange variant
927  element::dpc_variant _dpc_variant;
928 
929  // Degree that was input when creating the element
930  int _degree;
931 
932  // Highest degree polynomial in element's polyset
933  int _highest_degree;
934 
935  // Highest degree space that is a subspace of element's polyset
936  int _highest_complete_degree;
937 
938  // Value shape
939  std::vector<int> _value_shape;
940 
942  maps::type _map_type;
943 
944  // Shape function coefficient of expansion sets on cell. If shape
945  // function is given by @f$\psi_i = \sum_{k} \phi_{k}
946  // \alpha^{i}_{k}@f$, then _coeffs(i, j) = @f$\alpha^i_k@f$. ie
947  // _coeffs.row(i) are the expansion coefficients for shape function i
948  // (@f$\psi_{i}@f$).
949  xt::xtensor<double, 2> _coeffs;
950 
951  // Number of dofs associated with each cell (sub-)entity
952  //
953  // The dofs of an element are associated with entities of different
954  // topological dimension (vertices, edges, faces, cells). The dofs are
955  // listed in this order, with vertex dofs first. Each entry is the dof
956  // count on the associated entity, as listed by cell::topology.
957  std::vector<std::vector<int>> _num_edofs;
958 
959  // Number of dofs associated with the closure of each cell
960  // (sub-)entity
961  std::vector<std::vector<int>> _num_e_closure_dofs;
962 
963  // Dofs associated with each cell (sub-)entity
964  std::vector<std::vector<std::vector<int>>> _edofs;
965 
966  // Dofs associated with each cell (sub-)entity
967  std::vector<std::vector<std::vector<int>>> _e_closure_dofs;
968 
969  // Entity transformations
970  std::map<cell::type, xt::xtensor<double, 3>> _entity_transformations;
971 
972  // Set of points used for point evaluation
973  // Experimental - currently used for an implementation of
974  // "tabulate_dof_coordinates" Most useful for Lagrange. This may change or go
975  // away. For non-Lagrange elements, these points will be used in combination
976  // with _interpolation_matrix to perform interpolation
977  xt::xtensor<double, 2> _points;
978 
979  // Interpolation points on the cell. The shape is (entity_dim, num
980  // entities of given dimension, num_points, tdim)
981  std::array<std::vector<xt::xtensor<double, 2>>, 4> _x;
982 
984  xt::xtensor<double, 2> _matM;
985 
986  // Indicates whether or not the DOF transformations are all
987  // permutations
988  bool _dof_transformations_are_permutations;
989 
990  // Indicates whether or not the DOF transformations are all identity
991  bool _dof_transformations_are_identity;
992 
993  // The entity permutations (factorised). This will only be set if
994  // _dof_transformations_are_permutations is True and
995  // _dof_transformations_are_identity is False
996  std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm;
997 
998  // The reverse entity permutations (factorised). This will only be set
999  // if _dof_transformations_are_permutations is True and
1000  // _dof_transformations_are_identity is False
1001  std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm_rev;
1002 
1003  // The entity transformations in precomputed form
1004  std::map<cell::type,
1005  std::vector<std::tuple<std::vector<std::size_t>, std::vector<double>,
1006  xt::xtensor<double, 2>>>>
1007  _etrans;
1008 
1009  // The transposed entity transformations in precomputed form
1010  std::map<cell::type,
1011  std::vector<std::tuple<std::vector<std::size_t>, std::vector<double>,
1012  xt::xtensor<double, 2>>>>
1013  _etransT;
1014 
1015  // The inverse entity transformations in precomputed form
1016  std::map<cell::type,
1017  std::vector<std::tuple<std::vector<std::size_t>, std::vector<double>,
1018  xt::xtensor<double, 2>>>>
1019  _etrans_inv;
1020 
1021  // The inverse transpose entity transformations in precomputed form
1022  std::map<cell::type,
1023  std::vector<std::tuple<std::vector<std::size_t>, std::vector<double>,
1024  xt::xtensor<double, 2>>>>
1025  _etrans_invT;
1026 
1027  // Indicates whether or not this is the discontinuous version of the
1028  // element
1029  bool _discontinuous;
1030 
1031  // The dual matrix
1032  xt::xtensor<double, 2> _dual_matrix;
1033 
1034  // Tensor product representation
1035  // Entries of tuple are (list of elements on an interval, permutation
1036  // of DOF numbers)
1037  // @todo: For vector-valued elements, a tensor product type and a
1038  // scaling factor may additionally be needed.
1039  std::vector<std::tuple<std::vector<FiniteElement>, std::vector<int>>>
1040  _tensor_factors;
1041 
1042  // Is the interpolation matrix an identity?
1043  bool _interpolation_is_identity;
1044 
1045  // The coefficients that define the polynomial set in terms of the orthonormal
1046  // polynomials
1047  xt::xtensor<double, 2> _wcoeffs;
1048 
1049  // Interpolation matrices for each entity
1050  std::array<std::vector<xt::xtensor<double, 3>>, 4> _M;
1051 };
1052 
1073 FiniteElement create_custom_element(
1074  cell::type cell_type, const std::vector<std::size_t>& value_shape,
1075  const xt::xtensor<double, 2>& wcoeffs,
1076  const std::array<std::vector<xt::xtensor<double, 2>>, 4>& x,
1077  const std::array<std::vector<xt::xtensor<double, 3>>, 4>& M,
1078  maps::type map_type, bool discontinuous, int highest_complete_degree,
1079  int highest_degree);
1080 
1090 FiniteElement create_element(element::family family, cell::type cell,
1091  int degree, element::lagrange_variant lvariant,
1092  bool discontinuous);
1093 
1104 FiniteElement create_element(element::family family, cell::type cell,
1105  int degree, element::lagrange_variant lvariant,
1106  element::dpc_variant dvariant, bool discontinuous);
1107 
1117 FiniteElement create_element(element::family family, cell::type cell,
1118  int degree, element::dpc_variant dvariant,
1119  bool discontinuous);
1120 
1131 FiniteElement create_element(element::family family, cell::type cell,
1132  int degree, bool discontinuous);
1133 
1141 FiniteElement create_element(element::family family, cell::type cell,
1142  int degree, element::lagrange_variant lvariant);
1143 
1153 FiniteElement create_element(element::family family, cell::type cell,
1154  int degree, element::lagrange_variant lvariant,
1155  element::dpc_variant dvariant);
1156 
1164 FiniteElement create_element(element::family family, cell::type cell,
1165  int degree, element::dpc_variant dvariant);
1166 
1173 FiniteElement create_element(element::family family, cell::type cell,
1174  int degree);
1175 
1178 std::string version();
1179 
1180 //-----------------------------------------------------------------------------
1181 template <typename T>
1182 void FiniteElement::apply_dof_transformation(const xtl::span<T>& data,
1183  int block_size,
1184  std::uint32_t cell_info) const
1185 {
1186  if (_dof_transformations_are_identity)
1187  return;
1188 
1189  if (_cell_tdim >= 2)
1190  {
1191  // This assumes 3 bits are used per face. This will need updating if
1192  // 3D cells with faces with more than 4 sides are implemented
1193  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1194  int dofstart
1195  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1196 
1197  // Transform DOFs on edges
1198  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1199  {
1200  // Reverse an edge
1201  if (cell_info >> (face_start + e) & 1)
1202  precompute::apply_matrix(_etrans.at(cell::type::interval)[0], data,
1203  dofstart, block_size);
1204  dofstart += _num_edofs[1][e];
1205  }
1206 
1207  if (_cell_tdim == 3)
1208  {
1209  // Permute DOFs on faces
1210  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1211  {
1212  // Reflect a face
1213  if (cell_info >> (3 * f) & 1)
1214  precompute::apply_matrix(_etrans.at(_cell_subentity_types[2][f])[1],
1215  data, dofstart, block_size);
1216 
1217  // Rotate a face
1218  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1219  precompute::apply_matrix(_etrans.at(_cell_subentity_types[2][f])[0],
1220  data, dofstart, block_size);
1221  dofstart += _num_edofs[2][f];
1222  }
1223  }
1224  }
1225 }
1226 //-----------------------------------------------------------------------------
1227 template <typename T>
1229  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1230 {
1231  if (_dof_transformations_are_identity)
1232  return;
1233 
1234  if (_cell_tdim >= 2)
1235  {
1236  // This assumes 3 bits are used per face. This will need updating if
1237  // 3D cells with faces with more than 4 sides are implemented
1238  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1239  int dofstart
1240  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1241 
1242  // Transform DOFs on edges
1243  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1244  {
1245  // Reverse an edge
1246  if (cell_info >> (face_start + e) & 1)
1247  precompute::apply_matrix(_etransT.at(cell::type::interval)[0], data,
1248  dofstart, block_size);
1249  dofstart += _num_edofs[1][e];
1250  }
1251 
1252  if (_cell_tdim == 3)
1253  {
1254  // Permute DOFs on faces
1255  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1256  {
1257  // Rotate a face
1258  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1259  precompute::apply_matrix(_etransT.at(_cell_subentity_types[2][f])[0],
1260  data, dofstart, block_size);
1261  // Reflect a face
1262  if (cell_info >> (3 * f) & 1)
1263  precompute::apply_matrix(_etransT.at(_cell_subentity_types[2][f])[1],
1264  data, dofstart, block_size);
1265  dofstart += _num_edofs[2][f];
1266  }
1267  }
1268  }
1269 }
1270 //-----------------------------------------------------------------------------
1271 template <typename T>
1273  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1274 {
1275  if (_dof_transformations_are_identity)
1276  return;
1277 
1278  if (_cell_tdim >= 2)
1279  {
1280  // This assumes 3 bits are used per face. This will need updating if
1281  // 3D cells with faces with more than 4 sides are implemented
1282  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1283  int dofstart
1284  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1285 
1286  // Transform DOFs on edges
1287  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1288  {
1289  // Reverse an edge
1290  if (cell_info >> (face_start + e) & 1)
1291  precompute::apply_matrix(_etrans_invT.at(cell::type::interval)[0], data,
1292  dofstart, block_size);
1293  dofstart += _num_edofs[1][e];
1294  }
1295 
1296  if (_cell_tdim == 3)
1297  {
1298  // Permute DOFs on faces
1299  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1300  {
1301  // Reflect a face
1302  if (cell_info >> (3 * f) & 1)
1304  _etrans_invT.at(_cell_subentity_types[2][f])[1], data, dofstart,
1305  block_size);
1306 
1307  // Rotate a face
1308  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1310  _etrans_invT.at(_cell_subentity_types[2][f])[0], data, dofstart,
1311  block_size);
1312  dofstart += _num_edofs[2][f];
1313  }
1314  }
1315  }
1316 }
1317 //-----------------------------------------------------------------------------
1318 template <typename T>
1320  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1321 {
1322  if (_dof_transformations_are_identity)
1323  return;
1324 
1325  if (_cell_tdim >= 2)
1326  {
1327  // This assumes 3 bits are used per face. This will need updating if
1328  // 3D cells with faces with more than 4 sides are implemented
1329  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1330  int dofstart
1331  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1332 
1333  // Transform DOFs on edges
1334  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1335  {
1336  // Reverse an edge
1337  if (cell_info >> (face_start + e) & 1)
1338  precompute::apply_matrix(_etrans_inv.at(cell::type::interval)[0], data,
1339  dofstart, block_size);
1340  dofstart += _num_edofs[1][e];
1341  }
1342 
1343  if (_cell_tdim == 3)
1344  {
1345  // Permute DOFs on faces
1346  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1347  {
1348  // Rotate a face
1349  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1351  _etrans_inv.at(_cell_subentity_types[2][f])[0], data, dofstart,
1352  block_size);
1353  // Reflect a face
1354  if (cell_info >> (3 * f) & 1)
1356  _etrans_inv.at(_cell_subentity_types[2][f])[1], data, dofstart,
1357  block_size);
1358  dofstart += _num_edofs[2][f];
1359  }
1360  }
1361  }
1362 }
1363 //-----------------------------------------------------------------------------
1364 template <typename T>
1366  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1367 {
1368  if (_dof_transformations_are_identity)
1369  return;
1370 
1371  if (_cell_tdim >= 2)
1372  {
1373  // This assumes 3 bits are used per face. This will need updating if
1374  // 3D cells with faces with more than 4 sides are implemented
1375  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1376  int dofstart
1377  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1378 
1379  // Transform DOFs on edges
1380  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1381  {
1382  // Reverse an edge
1383  if (cell_info >> (face_start + e) & 1)
1385  _etrans.at(cell::type::interval)[0], data, dofstart, block_size);
1386  dofstart += _num_edofs[1][e];
1387  }
1388 
1389  if (_cell_tdim == 3)
1390  {
1391  // Permute DOFs on faces
1392  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1393  {
1394  // Reflect a face
1395  if (cell_info >> (3 * f) & 1)
1397  _etrans.at(_cell_subentity_types[2][f])[1], data, dofstart,
1398  block_size);
1399 
1400  // Rotate a face
1401  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1403  _etrans.at(_cell_subentity_types[2][f])[0], data, dofstart,
1404  block_size);
1405  dofstart += _num_edofs[2][f];
1406  }
1407  }
1408  }
1409 }
1410 //-----------------------------------------------------------------------------
1411 template <typename T>
1413  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1414 {
1415  if (_dof_transformations_are_identity)
1416  return;
1417 
1418  if (_cell_tdim >= 2)
1419  {
1420  // This assumes 3 bits are used per face. This will need updating if
1421  // 3D cells with faces with more than 4 sides are implemented
1422  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1423  int dofstart
1424  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1425 
1426  // Transform DOFs on edges
1427  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1428  {
1429  // Reverse an edge
1430  if (cell_info >> (face_start + e) & 1)
1432  _etrans_invT.at(cell::type::interval)[0], data, dofstart,
1433  block_size);
1434  dofstart += _num_edofs[1][e];
1435  }
1436 
1437  if (_cell_tdim == 3)
1438  {
1439  // Permute DOFs on faces
1440  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1441  {
1442  // Reflect a face
1443  if (cell_info >> (3 * f) & 1)
1445  _etrans_invT.at(_cell_subentity_types[2][f])[1], data, dofstart,
1446  block_size);
1447 
1448  // Rotate a face
1449  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1451  _etrans_invT.at(_cell_subentity_types[2][f])[0], data, dofstart,
1452  block_size);
1453  dofstart += _num_edofs[2][f];
1454  }
1455  }
1456  }
1457 }
1458 //-----------------------------------------------------------------------------
1459 template <typename T>
1461  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1462 {
1463  if (_dof_transformations_are_identity)
1464  return;
1465 
1466  if (_cell_tdim >= 2)
1467  {
1468  // This assumes 3 bits are used per face. This will need updating if
1469  // 3D cells with faces with more than 4 sides are implemented
1470  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1471  int dofstart
1472  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1473 
1474  // Transform DOFs on edges
1475  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1476  {
1477  // Reverse an edge
1478  if (cell_info >> (face_start + e) & 1)
1480  _etransT.at(cell::type::interval)[0], data, dofstart, block_size);
1481  dofstart += _num_edofs[1][e];
1482  }
1483 
1484  if (_cell_tdim == 3)
1485  {
1486  // Permute DOFs on faces
1487  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1488  {
1489  // Rotate a face
1490  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1492  _etransT.at(_cell_subentity_types[2][f])[0], data, dofstart,
1493  block_size);
1494 
1495  // Reflect a face
1496  if (cell_info >> (3 * f) & 1)
1498  _etransT.at(_cell_subentity_types[2][f])[1], data, dofstart,
1499  block_size);
1500 
1501  dofstart += _num_edofs[2][f];
1502  }
1503  }
1504  }
1505 }
1506 //-----------------------------------------------------------------------------
1507 template <typename T>
1509  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1510 {
1511  if (_dof_transformations_are_identity)
1512  return;
1513 
1514  if (_cell_tdim >= 2)
1515  {
1516  // This assumes 3 bits are used per face. This will need updating if
1517  // 3D cells with faces with more than 4 sides are implemented
1518  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1519  int dofstart
1520  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1521 
1522  // Transform DOFs on edges
1523  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1524  {
1525  // Reverse an edge
1526  if (cell_info >> (face_start + e) & 1)
1528  _etrans_inv.at(cell::type::interval)[0], data, dofstart,
1529  block_size);
1530  dofstart += _num_edofs[1][e];
1531  }
1532 
1533  if (_cell_tdim == 3)
1534  {
1535  // Permute DOFs on faces
1536  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1537  {
1538  // Rotate a face
1539  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1541  _etrans_inv.at(_cell_subentity_types[2][f])[0], data, dofstart,
1542  block_size);
1543 
1544  // Reflect a face
1545  if (cell_info >> (3 * f) & 1)
1547  _etrans_inv.at(_cell_subentity_types[2][f])[1], data, dofstart,
1548  block_size);
1549 
1550  dofstart += _num_edofs[2][f];
1551  }
1552  }
1553  }
1554 }
1555 //-----------------------------------------------------------------------------
1556 
1557 } // namespace basix
basix::FiniteElement::map_type
maps::type map_type() const
Definition: finite-element.cpp:877
basix::FiniteElement::permute_dofs
void permute_dofs(const xtl::span< std::int32_t > &dofs, std::uint32_t cell_info) const
Definition: finite-element.cpp:1035
basix::FiniteElement::highest_complete_degree
int highest_complete_degree() const
Definition: finite-element.cpp:863
basix::FiniteElement::dof_transformations_are_identity
bool dof_transformations_are_identity() const
Definition: finite-element.cpp:886
basix::FiniteElement::M
const std::array< std::vector< xt::xtensor< double, 3 > >, 4 > & M() const
Definition: finite-element.cpp:1171
basix::FiniteElement::lagrange_variant
element::lagrange_variant lagrange_variant() const
Definition: finite-element.cpp:1188
basix::element::lagrange_variant
lagrange_variant
Variants of a Lagrange space that can be created.
Definition: element-families.h:11
basix::FiniteElement::apply_transpose_dof_transformation
void apply_transpose_dof_transformation(const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
Definition: finite-element.h:1228
basix::FiniteElement::apply_dof_transformation
void apply_dof_transformation(const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
Definition: finite-element.h:1182
basix::FiniteElement::degree
int degree() const
Definition: finite-element.cpp:859
basix::FiniteElement::~FiniteElement
~FiniteElement()=default
Destructor.
basix::FiniteElement::dpc_variant
element::dpc_variant dpc_variant() const
Definition: finite-element.cpp:1193
basix::maps::l2_piola
void l2_piola(O &&r, const P &U, const Q &, double detJ, const R &)
L2 Piola map.
Definition: maps.h:63
basix::FiniteElement::num_entity_dofs
const std::vector< std::vector< int > > & num_entity_dofs() const
Definition: finite-element.cpp:896
basix::cell::type
type
Cell type.
Definition: cell.h:18
basix::FiniteElement::points
const xt::xtensor< double, 2 > & points() const
Definition: finite-element.cpp:976
basix::FiniteElement::push_forward
xt::xtensor< double, 3 > push_forward(const xt::xtensor< double, 3 > &U, const xt::xtensor< double, 3 > &J, const xtl::span< const double > &detJ, const xt::xtensor< double, 3 > &K) const
Definition: finite-element.cpp:978
basix::FiniteElement::has_tensor_product_factorisation
bool has_tensor_product_factorisation() const
Definition: finite-element.cpp:1183
basix::FiniteElement::tabulate
xt::xtensor< double, 4 > tabulate(int nd, const xt::xtensor< double, 2 > &x) const
Definition: finite-element.cpp:819
basix::FiniteElement::entity_transformations
std::map< cell::type, xt::xtensor< double, 3 > > entity_transformations() const
Definition: finite-element.cpp:1147
basix::FiniteElement::cell_type
cell::type cell_type() const
Definition: finite-element.cpp:857
basix::FiniteElement::apply_inverse_transpose_dof_transformation_to_transpose
void apply_inverse_transpose_dof_transformation_to_transpose(const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
Definition: finite-element.h:1412
basix::FiniteElement::entity_dofs
const std::vector< std::vector< std::vector< int > > > & entity_dofs() const
Definition: finite-element.cpp:902
basix::maps::double_contravariant_piola
void double_contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Double contravariant Piola map.
Definition: maps.h:115
basix::FiniteElement::discontinuous
bool discontinuous() const
Definition: finite-element.cpp:879
basix::FiniteElement::apply_inverse_dof_transformation
void apply_inverse_dof_transformation(const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
Definition: finite-element.h:1319
basix::FiniteElement::apply_transpose_dof_transformation_to_transpose
void apply_transpose_dof_transformation_to_transpose(const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
Definition: finite-element.h:1460
basix::create_element
FiniteElement create_element(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, bool discontinuous)
Definition: finite-element.cpp:289
basix::FiniteElement::dof_transformations_are_permutations
bool dof_transformations_are_permutations() const
Definition: finite-element.cpp:881
basix::FiniteElement::family
element::family family() const
Definition: finite-element.cpp:875
basix::FiniteElement::apply_inverse_dof_transformation_to_transpose
void apply_inverse_dof_transformation_to_transpose(const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
Definition: finite-element.h:1508
basix::FiniteElement::get_tensor_product_representation
std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int > > > get_tensor_product_representation() const
Definition: finite-element.cpp:1201
basix::element::make_discontinuous
std::tuple< std::array< std::vector< xt::xtensor< double, 2 > >, 4 >, std::array< std::vector< xt::xtensor< double, 3 > >, 4 > > make_discontinuous(const std::array< std::vector< xt::xtensor< double, 2 >>, 4 > &x, const std::array< std::vector< xt::xtensor< double, 3 >>, 4 > &M, int tdim, int value_size)
Definition: finite-element.cpp:336
basix::FiniteElement::dim
int dim() const
Definition: finite-element.cpp:873
basix::FiniteElement::operator==
bool operator==(const FiniteElement &e) const
Definition: finite-element.cpp:785
basix::FiniteElement::apply_dof_transformation_to_transpose
void apply_dof_transformation_to_transpose(const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
Definition: finite-element.h:1365
basix::FiniteElement::tabulate_shape
std::array< std::size_t, 4 > tabulate_shape(std::size_t nd, std::size_t num_points) const
Definition: finite-element.cpp:805
basix::FiniteElement::highest_degree
int highest_degree() const
Definition: finite-element.cpp:861
basix::maps::type
type
Map type.
Definition: maps.h:19
basix::FiniteElement::value_shape
const std::vector< int > & value_shape() const
Definition: finite-element.cpp:868
basix::FiniteElement::base_transformations
xt::xtensor< double, 3 > base_transformations() const
Definition: finite-element.cpp:919
basix::FiniteElement::interpolation_is_identity
bool interpolation_is_identity() const
Definition: finite-element.cpp:1195
basix::FiniteElement::num_entity_closure_dofs
const std::vector< std::vector< int > > & num_entity_closure_dofs() const
Definition: finite-element.cpp:908
basix::maps::double_covariant_piola
void double_covariant_piola(O &&r, const P &U, const Q &J, double, const R &K)
Double covariant Piola map.
Definition: maps.h:100
basix::precompute::apply_matrix_to_transpose
void apply_matrix_to_transpose(const std::tuple< std::vector< std::size_t >, std::vector< T >, xt::xtensor< T, 2 >> &matrix, const xtl::span< E > &data, std::size_t offset=0, std::size_t block_size=1)
Definition: precompute.h:351
basix::FiniteElement
A finite element.
Definition: finite-element.h:52
basix::FiniteElement::interpolation_matrix
const xt::xtensor< double, 2 > & interpolation_matrix() const
Definition: finite-element.cpp:891
basix::FiniteElement::wcoeffs
const xt::xtensor< double, 2 > & wcoeffs() const
Definition: finite-element.cpp:1157
basix::maps::contravariant_piola
void contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Contravariant Piola map.
Definition: maps.h:85
basix::FiniteElement::FiniteElement
FiniteElement(element::family family, cell::type cell_type, int degree, const std::vector< std::size_t > &value_shape, const xt::xtensor< double, 2 > &wcoeffs, const std::array< std::vector< xt::xtensor< double, 2 >>, 4 > &x, const std::array< std::vector< xt::xtensor< double, 3 >>, 4 > &M, maps::type map_type, bool discontinuous, int highest_complete_degree, int highest_degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, std::vector< std::tuple< std::vector< FiniteElement >, std::vector< int >>> tensor_factors={})
Definition: finite-element.cpp:516
basix::FiniteElement::coefficient_matrix
const xt::xtensor< double, 2 > & coefficient_matrix() const
Definition: finite-element.cpp:1178
basix::FiniteElement::dual_matrix
const xt::xtensor< double, 2 > & dual_matrix() const
Definition: finite-element.cpp:1152
basix::version
std::string version()
Definition: finite-element.cpp:1208
basix::FiniteElement::entity_closure_dofs
const std::vector< std::vector< std::vector< int > > > & entity_closure_dofs() const
Definition: finite-element.cpp:914
basix::FiniteElement::map_fn
std::function< void(O &, const P &, const Q &, double, const R &)> map_fn() const
Definition: finite-element.h:477
basix::create_custom_element
FiniteElement create_custom_element(cell::type cell_type, const std::vector< std::size_t > &value_shape, const xt::xtensor< double, 2 > &wcoeffs, const std::array< std::vector< xt::xtensor< double, 2 >>, 4 > &x, const std::array< std::vector< xt::xtensor< double, 3 >>, 4 > &M, maps::type map_type, bool discontinuous, int highest_complete_degree, int highest_degree)
Definition: finite-element.cpp:392
basix
Basix: FEniCS runtime basis evaluation library.
Definition: cell.h:14
basix::FiniteElement::unpermute_dofs
void unpermute_dofs(const xtl::span< std::int32_t > &dofs, std::uint32_t cell_info) const
Definition: finite-element.cpp:1091
basix::precompute::apply_matrix
void apply_matrix(const std::tuple< std::vector< std::size_t >, std::vector< T >, xt::xtensor< T, 2 >> &matrix, const xtl::span< E > &data, std::size_t offset=0, std::size_t block_size=1)
Definition: precompute.h:319
basix::FiniteElement::operator=
FiniteElement & operator=(const FiniteElement &element)=default
Assignment operator.
basix::FiniteElement::x
const std::array< std::vector< xt::xtensor< double, 2 > >, 4 > & x() const
Definition: finite-element.cpp:1165
basix::FiniteElement::pull_back
xt::xtensor< double, 3 > pull_back(const xt::xtensor< double, 3 > &u, const xt::xtensor< double, 3 > &J, const xtl::span< const double > &detJ, const xt::xtensor< double, 3 > &K) const
Definition: finite-element.cpp:1006
basix::maps::covariant_piola
void covariant_piola(O &&r, const P &U, const Q &, double, const R &K)
Covariant Piola map.
Definition: maps.h:71
basix::element::dpc_variant
dpc_variant
Definition: element-families.h:31
basix::FiniteElement::apply_inverse_transpose_dof_transformation
void apply_inverse_transpose_dof_transformation(const xtl::span< T > &data, int block_size, std::uint32_t cell_info) const
Definition: finite-element.h:1272
basix::element::family
family
Available element families.
Definition: element-families.h:44