Basix 0.8.0

Home     Installation     Demos     C++ docs     Python docs

Loading...
Searching...
No Matches
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
27namespace basix
28{
29
30namespace impl
31{
32template <typename T, std::size_t d>
33using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
34 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
35template <typename T, std::size_t d>
36using mdarray_t
37 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::mdarray<
38 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
39
42template <typename T>
43std::array<std::vector<mdspan_t<const T, 2>>, 4>
44to_mdspan(std::array<std::vector<mdarray_t<T, 2>>, 4>& x)
45{
46 std::array<std::vector<mdspan_t<const T, 2>>, 4> x1;
47 for (std::size_t i = 0; i < x.size(); ++i)
48 for (std::size_t j = 0; j < x[i].size(); ++j)
49 x1[i].emplace_back(x[i][j].data(), x[i][j].extents());
50
51 return x1;
52}
53
56template <typename T>
57std::array<std::vector<mdspan_t<const T, 4>>, 4>
58to_mdspan(std::array<std::vector<mdarray_t<T, 4>>, 4>& M)
59{
60 std::array<std::vector<mdspan_t<const T, 4>>, 4> M1;
61 for (std::size_t i = 0; i < M.size(); ++i)
62 for (std::size_t j = 0; j < M[i].size(); ++j)
63 M1[i].emplace_back(M[i][j].data(), M[i][j].extents());
64
65 return M1;
66}
67
70template <typename T>
71std::array<std::vector<mdspan_t<const T, 2>>, 4>
72to_mdspan(const std::array<std::vector<std::vector<T>>, 4>& x,
73 const std::array<std::vector<std::array<std::size_t, 2>>, 4>& shape)
74{
75 std::array<std::vector<mdspan_t<const T, 2>>, 4> x1;
76 for (std::size_t i = 0; i < x.size(); ++i)
77 for (std::size_t j = 0; j < x[i].size(); ++j)
78 x1[i].push_back(mdspan_t<const T, 2>(x[i][j].data(), shape[i][j]));
79
80 return x1;
81}
82
85template <typename T>
86std::array<std::vector<mdspan_t<const T, 4>>, 4>
87to_mdspan(const std::array<std::vector<std::vector<T>>, 4>& M,
88 const std::array<std::vector<std::array<std::size_t, 4>>, 4>& shape)
89{
90 std::array<std::vector<mdspan_t<const T, 4>>, 4> M1;
91 for (std::size_t i = 0; i < M.size(); ++i)
92 for (std::size_t j = 0; j < M[i].size(); ++j)
93 M1[i].push_back(mdspan_t<const T, 4>(M[i][j].data(), shape[i][j]));
94
95 return M1;
96}
97
98} // namespace impl
99
100namespace element
101{
103template <typename T, std::size_t d>
104using mdspan_t = impl::mdspan_t<T, d>;
105
121template <std::floating_point T>
122std::tuple<std::array<std::vector<std::vector<T>>, 4>,
123 std::array<std::vector<std::array<std::size_t, 2>>, 4>,
124 std::array<std::vector<std::vector<T>>, 4>,
125 std::array<std::vector<std::array<std::size_t, 4>>, 4>>
126make_discontinuous(const std::array<std::vector<mdspan_t<const T, 2>>, 4>& x,
127 const std::array<std::vector<mdspan_t<const T, 4>>, 4>& M,
128 std::size_t tdim, std::size_t value_size);
129
130} // namespace element
131
137template <std::floating_point F>
139{
140 template <typename T, std::size_t d>
141 using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
142 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, d>>;
143
144public:
318 polyset::type poly_type, int degree,
319 const std::vector<std::size_t>& value_shape,
320 mdspan_t<const F, 2> wcoeffs,
321 const std::array<std::vector<mdspan_t<const F, 2>>, 4>& x,
322 const std::array<std::vector<mdspan_t<const F, 4>>, 4>& M,
327 element::dpc_variant dvariant,
328 std::vector<int> dof_ordering = {});
329
331 FiniteElement(const FiniteElement& element) = default;
332
334 FiniteElement(FiniteElement&& element) = default;
335
337 ~FiniteElement() = default;
338
340 FiniteElement& operator=(const FiniteElement& element) = default;
341
343 FiniteElement& operator=(FiniteElement&& element) = default;
344
349 bool operator==(const FiniteElement& e) const;
350
352 std::size_t hash() const;
353
363 std::array<std::size_t, 4> tabulate_shape(std::size_t nd,
364 std::size_t num_points) const
365 {
366 std::size_t ndsize = 1;
367 for (std::size_t i = 1; i <= nd; ++i)
368 ndsize *= (_cell_tdim + i);
369 for (std::size_t i = 1; i <= nd; ++i)
370 ndsize /= i;
371 std::size_t vs = std::accumulate(_value_shape.begin(), _value_shape.end(),
372 1, std::multiplies{});
373 std::size_t ndofs = _coeffs.second[0];
374 return {ndsize, num_points, ndofs, vs};
375 }
376
398 std::pair<std::vector<F>, std::array<std::size_t, 4>>
399 tabulate(int nd, impl::mdspan_t<const F, 2> x) const;
400
424 std::pair<std::vector<F>, std::array<std::size_t, 4>>
425 tabulate(int nd, std::span<const F> x,
426 std::array<std::size_t, 2> shape) const;
427
453 void tabulate(int nd, impl::mdspan_t<const F, 2> x,
454 mdspan_t<F, 4> basis) const;
455
481 void tabulate(int nd, std::span<const F> x, std::array<std::size_t, 2> xshape,
482 std::span<F> basis) const;
483
486 cell::type cell_type() const { return _cell_type; }
487
490 polyset::type polyset_type() const { return _poly_type; }
491
494 int degree() const { return _degree; }
495
500 int embedded_superdegree() const { return _embedded_superdegree; }
501
505 int embedded_subdegree() const { return _embedded_subdegree; }
506
512 const std::vector<std::size_t>& value_shape() const { return _value_shape; }
513
518 int dim() const { return _coeffs.second[0]; }
519
522 element::family family() const { return _family; }
523
527 {
528 return _lagrange_variant;
529 }
530
533 element::dpc_variant dpc_variant() const { return _dpc_variant; }
534
537 maps::type map_type() const { return _map_type; }
538
541 sobolev::space sobolev_space() const { return _sobolev_space; }
542
547 bool discontinuous() const { return _discontinuous; }
548
552 {
553 return _dof_transformations_are_permutations;
554 }
555
558 {
559 return _dof_transformations_are_identity;
560 }
561
576 std::pair<std::vector<F>, std::array<std::size_t, 3>>
577 push_forward(impl::mdspan_t<const F, 3> U, impl::mdspan_t<const F, 3> J,
578 std::span<const F> detJ, impl::mdspan_t<const F, 3> K) const;
579
587 std::pair<std::vector<F>, std::array<std::size_t, 3>>
588 pull_back(impl::mdspan_t<const F, 3> u, impl::mdspan_t<const F, 3> J,
589 std::span<const F> detJ, impl::mdspan_t<const F, 3> K) const;
590
622 template <typename O, typename P, typename Q, typename R>
623 std::function<void(O&, const P&, const Q&, F, const R&)> map_fn() const
624 {
625 switch (_map_type)
626 {
627 case maps::type::identity:
628 return [](O& u, const P& U, const Q&, F, const R&)
629 {
630 assert(U.extent(0) == u.extent(0));
631 assert(U.extent(1) == u.extent(1));
632 for (std::size_t i = 0; i < U.extent(0); ++i)
633 for (std::size_t j = 0; j < U.extent(1); ++j)
634 u(i, j) = U(i, j);
635 };
636 case maps::type::covariantPiola:
637 return [](O& u, const P& U, const Q& J, F detJ, const R& K)
638 { maps::covariant_piola(u, U, J, detJ, K); };
639 case maps::type::contravariantPiola:
640 return [](O& u, const P& U, const Q& J, F detJ, const R& K)
641 { maps::contravariant_piola(u, U, J, detJ, K); };
642 case maps::type::doubleCovariantPiola:
643 return [](O& u, const P& U, const Q& J, F detJ, const R& K)
644 { maps::double_covariant_piola(u, U, J, detJ, K); };
645 case maps::type::doubleContravariantPiola:
646 return [](O& u, const P& U, const Q& J, F detJ, const R& K)
647 { maps::double_contravariant_piola(u, U, J, detJ, K); };
648 default:
649 throw std::runtime_error("Map not implemented");
650 }
651 }
652
660 const std::vector<std::vector<std::vector<int>>>& entity_dofs() const
661 {
662 return _edofs;
663 }
664
674 const std::vector<std::vector<std::vector<int>>>& entity_closure_dofs() const
675 {
676 return _e_closure_dofs;
677 }
678
760 std::pair<std::vector<F>, std::array<std::size_t, 3>>
761 base_transformations() const;
762
767 std::map<cell::type, std::pair<std::vector<F>, std::array<std::size_t, 3>>>
769 {
770 return _entity_transformations;
771 }
772
795 void permute(std::span<std::int32_t> d, std::uint32_t cell_info) const
796 {
797 if (!_dof_transformations_are_permutations)
798 {
799 throw std::runtime_error(
800 "The DOF transformations for this element are not permutations");
801 }
802
803 if (_dof_transformations_are_identity)
804 return;
805 else
806 permute_data<std::int32_t, false>(d, 1, cell_info, _eperm);
807 }
808
826 void permute_inv(std::span<std::int32_t> d, std::uint32_t cell_info) const
827 {
828 if (!_dof_transformations_are_permutations)
829 {
830 throw std::runtime_error(
831 "The DOF transformations for this element are not permutations");
832 }
833
834 if (_dof_transformations_are_identity)
835 return;
836 else
837 permute_data<std::int32_t, true>(d, 1, cell_info, _eperm_rev);
838 }
839
870 template <typename T>
871 void T_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
872
883 template <typename T>
884 void Tt_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
885
897 template <typename T>
898 void Tt_inv_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
899
910 template <typename T>
911 void Tinv_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
912
923 template <typename T>
924 void Tt_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
925
935 template <typename T>
936 void T_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
937
948 template <typename T>
949 void Tinv_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
950
961 template <typename T>
962 void Tt_inv_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
963
970 const std::pair<std::vector<F>, std::array<std::size_t, 2>>& points() const
971 {
972 return _points;
973 }
974
1026 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1028 {
1029 return _matM;
1030 }
1031
1037 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1039 {
1040 return _dual_matrix;
1041 }
1042
1079 const std::pair<std::vector<F>, std::array<std::size_t, 2>>& wcoeffs() const
1080 {
1081 return _wcoeffs;
1082 }
1083
1088 const std::array<
1089 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>, 4>&
1090 x() const
1091 {
1092 return _x;
1093 }
1094
1131 const std::array<
1132 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>, 4>&
1133 M() const
1134 {
1135 return _M;
1136 }
1137
1141 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1143 {
1144 return _coeffs;
1145 }
1146
1158 {
1159 return !_tensor_factors.empty();
1160 }
1161
1174 std::vector<std::vector<FiniteElement<F>>>
1176 {
1178 throw std::runtime_error("Element has no tensor product representation.");
1179 return _tensor_factors;
1180 }
1181
1186 bool interpolation_is_identity() const { return _interpolation_is_identity; }
1187
1189 int interpolation_nderivs() const { return _interpolation_nderivs; }
1190
1192 const std::vector<int>& dof_ordering() const { return _dof_ordering; }
1193
1194private:
1201 template <typename T, bool post>
1202 void permute_data(
1203 std::span<T> data, int block_size, std::uint32_t cell_info,
1204 const std::map<cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1205 const;
1206
1207 using array2_t = std::pair<std::vector<F>, std::array<std::size_t, 2>>;
1208 using array3_t = std::pair<std::vector<F>, std::array<std::size_t, 3>>;
1209 using trans_data_t
1210 = std::vector<std::pair<std::vector<std::size_t>, array2_t>>;
1211
1218 template <typename T, bool post, typename OP>
1219 void
1220 transform_data(std::span<T> data, int block_size, std::uint32_t cell_info,
1221 const std::map<cell::type, trans_data_t>& etrans, OP op) const;
1222
1223 // Cell type
1224 cell::type _cell_type;
1225
1226 // Polyset type
1227 polyset::type _poly_type;
1228
1229 // Topological dimension of the cell
1230 std::size_t _cell_tdim;
1231
1232 // Topological dimension of the cell
1233 std::vector<std::vector<cell::type>> _cell_subentity_types;
1234
1235 // Finite element family
1236 element::family _family;
1237
1238 // Lagrange variant
1239 element::lagrange_variant _lagrange_variant;
1240
1241 // DPC variant
1242 element::dpc_variant _dpc_variant;
1243
1244 // Degree that was input when creating the element
1245 int _degree;
1246
1247 // Degree
1248 int _interpolation_nderivs;
1249
1250 // Highest degree polynomial in element's polyset
1251 int _embedded_superdegree;
1252
1253 // Highest degree space that is a subspace of element's polyset
1254 int _embedded_subdegree;
1255
1256 // Value shape
1257 std::vector<std::size_t> _value_shape;
1258
1260 maps::type _map_type;
1261
1263 sobolev::space _sobolev_space;
1264
1265 // Shape function coefficient of expansion sets on cell. If shape
1266 // function is given by @f$\psi_i = \sum_{k} \phi_{k}
1267 // \alpha^{i}_{k}@f$, then _coeffs(i, j) = @f$\alpha^i_k@f$. ie
1268 // _coeffs.row(i) are the expansion coefficients for shape function i
1269 // (@f$\psi_{i}@f$).
1270 std::pair<std::vector<F>, std::array<std::size_t, 2>> _coeffs;
1271
1272 // Dofs associated with each cell (sub-)entity
1273 std::vector<std::vector<std::vector<int>>> _edofs;
1274
1275 // Dofs associated with the closdure of each cell (sub-)entity
1276 std::vector<std::vector<std::vector<int>>> _e_closure_dofs;
1277
1278 // Entity transformations
1279 std::map<cell::type, array3_t> _entity_transformations;
1280
1281 // Set of points used for point evaluation
1282 // Experimental - currently used for an implementation of
1283 // "tabulate_dof_coordinates" Most useful for Lagrange. This may change or go
1284 // away. For non-Lagrange elements, these points will be used in combination
1285 // with _interpolation_matrix to perform interpolation
1286 std::pair<std::vector<F>, std::array<std::size_t, 2>> _points;
1287
1288 // Interpolation points on the cell. The shape is (entity_dim, num
1289 // entities of given dimension, num_points, tdim)
1290 std::array<std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>,
1291 4>
1292 _x;
1293
1295 std::pair<std::vector<F>, std::array<std::size_t, 2>> _matM;
1296
1297 // Indicates whether or not the DOF transformations are all
1298 // permutations
1299 bool _dof_transformations_are_permutations;
1300
1301 // Indicates whether or not the DOF transformations are all identity
1302 bool _dof_transformations_are_identity;
1303
1304 // The entity permutations (factorised). This will only be set if
1305 // _dof_transformations_are_permutations is True and
1306 // _dof_transformations_are_identity is False
1307 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm;
1308
1309 // The reverse entity permutations (factorised). This will only be set
1310 // if _dof_transformations_are_permutations is True and
1311 // _dof_transformations_are_identity is False
1312 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm_rev;
1313
1314 // The entity transformations in precomputed form
1315 std::map<cell::type, trans_data_t> _etrans;
1316
1317 // The transposed entity transformations in precomputed form
1318 std::map<cell::type, trans_data_t> _etransT;
1319
1320 // The inverse entity transformations in precomputed form
1321 std::map<cell::type, trans_data_t> _etrans_inv;
1322
1323 // The inverse transpose entity transformations in precomputed form
1324 std::map<cell::type, trans_data_t> _etrans_invT;
1325
1326 // Indicates whether or not this is the discontinuous version of the
1327 // element
1328 bool _discontinuous;
1329
1330 // The dual matrix
1331 std::pair<std::vector<F>, std::array<std::size_t, 2>> _dual_matrix;
1332
1333 // Dof reordering for different element dof layout compatibility.
1334 // The reference basix layout is ordered by entity, i.e. dofs on
1335 // vertices, followed by edges, faces, then internal dofs.
1336 // _dof_ordering stores the map to the new order required, e.g.
1337 // for a P2 triangle, _dof_ordering=[0 3 5 1 2 4] will place
1338 // dofs 0, 3, 5 on the vertices and 1, 2, 4, on the edges.
1339 std::vector<int> _dof_ordering;
1340
1341 // Tensor product representation
1342 // Entries of tuple are (list of elements on an interval, permutation
1343 // of DOF numbers)
1344 // @todo: For vector-valued elements, a tensor product type and a
1345 // scaling factor may additionally be needed.
1346 std::vector<std::vector<FiniteElement>> _tensor_factors;
1347
1348 // Is the interpolation matrix an identity?
1349 bool _interpolation_is_identity;
1350
1351 // The coefficients that define the polynomial set in terms of the
1352 // orthonormal polynomials
1353 std::pair<std::vector<F>, std::array<std::size_t, 2>> _wcoeffs;
1354
1355 // Interpolation matrices for each entity
1356 using array4_t
1357 = std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>;
1358 std::array<array4_t, 4> _M;
1359};
1360
1386template <std::floating_point T>
1387FiniteElement<T> create_custom_element(
1388 cell::type cell_type, const std::vector<std::size_t>& value_shape,
1389 impl::mdspan_t<const T, 2> wcoeffs,
1390 const std::array<std::vector<impl::mdspan_t<const T, 2>>, 4>& x,
1391 const std::array<std::vector<impl::mdspan_t<const T, 4>>, 4>& M,
1392 int interpolation_nderivs, maps::type map_type,
1393 sobolev::space sobolev_space, bool discontinuous, int embedded_subdegree,
1394 int embedded_superdegree, polyset::type poly_type);
1395
1407template <std::floating_point T>
1408FiniteElement<T> create_element(element::family family, cell::type cell,
1409 int degree, element::lagrange_variant lvariant,
1410 element::dpc_variant dvariant,
1411 bool discontinuous,
1412 std::vector<int> dof_ordering = {});
1413
1425std::vector<int> tp_dof_ordering(element::family family, cell::type cell,
1426 int degree, element::lagrange_variant lvariant,
1427 element::dpc_variant dvariant,
1428 bool discontinuous);
1429
1442template <std::floating_point T>
1443std::vector<std::vector<FiniteElement<T>>>
1444tp_factors(element::family family, cell::type cell, int degree,
1446 bool discontinuous, std::vector<int> dof_ordering);
1447
1459template <std::floating_point T>
1460FiniteElement<T>
1461create_tp_element(element::family family, cell::type cell, int degree,
1463 element::dpc_variant dvariant, bool discontinuous);
1464
1467std::string version();
1468
1469//-----------------------------------------------------------------------------
1470template <std::floating_point F>
1471template <typename T, bool post>
1472void FiniteElement<F>::permute_data(
1473 std::span<T> data, int block_size, std::uint32_t cell_info,
1474 const std::map<cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1475 const
1476{
1477 if (_cell_tdim >= 2)
1478 {
1479 // This assumes 3 bits are used per face. This will need updating if 3D
1480 // cells with faces with more than 4 sides are implemented
1481 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1482
1483 // Permute DOFs on edges
1484 {
1485 auto& trans = eperm.at(cell::type::interval)[0];
1486 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1487 {
1488 // Reverse an edge
1489 if (cell_info >> (face_start + e) & 1)
1490 {
1491 precompute::apply_permutation_mapped(trans, data, _edofs[1][e],
1492 block_size);
1493 }
1494 }
1495 }
1496
1497 if (_cell_tdim == 3)
1498 {
1499 // Permute DOFs on faces
1500 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1501 {
1502 auto& trans = eperm.at(_cell_subentity_types[2][f]);
1503
1504 // Reflect a face (pre rotate)
1505 if (!post and cell_info >> (3 * f) & 1)
1506 {
1507 precompute::apply_permutation_mapped(trans[1], data, _edofs[2][f],
1508 block_size);
1509 }
1510
1511 // Rotate a face
1512 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1513 {
1514 precompute::apply_permutation_mapped(trans[0], data, _edofs[2][f],
1515 block_size);
1516 }
1517
1518 // Reflect a face (post rotate)
1519 if (post and cell_info >> (3 * f) & 1)
1520 {
1521 precompute::apply_permutation_mapped(trans[1], data, _edofs[2][f],
1522 block_size);
1523 }
1524 }
1525 }
1526 }
1527}
1528//-----------------------------------------------------------------------------
1529template <std::floating_point F>
1530template <typename T, bool post, typename OP>
1531void FiniteElement<F>::transform_data(
1532 std::span<T> data, int block_size, std::uint32_t cell_info,
1533 const std::map<cell::type, trans_data_t>& etrans, OP op) const
1534{
1535 if (_cell_tdim >= 2)
1536 {
1537 // This assumes 3 bits are used per face. This will need updating if
1538 // 3D cells with faces with more than 4 sides are implemented
1539 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1540 int dofstart = 0;
1541 for (auto& edofs0 : _edofs[0])
1542 dofstart += edofs0.size();
1543
1544 // Transform DOFs on edges
1545 {
1546 auto& [v_size_t, matrix] = etrans.at(cell::type::interval)[0];
1547 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1548 {
1549 // Reverse an edge
1550 if (cell_info >> (face_start + e) & 1)
1551 {
1552 op(std::span(v_size_t),
1553 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1554 dofstart, block_size);
1555 }
1556 dofstart += _edofs[1][e].size();
1557 }
1558 }
1559
1560 if (_cell_tdim == 3)
1561 {
1562 // Permute DOFs on faces
1563 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1564 {
1565 auto& trans = etrans.at(_cell_subentity_types[2][f]);
1566
1567 // Reflect a face (pre rotation)
1568 if (!post and cell_info >> (3 * f) & 1)
1569 {
1570 const auto& m = trans[1];
1571 const auto& v_size_t = std::get<0>(m);
1572 const auto& matrix = std::get<1>(m);
1573 op(std::span(v_size_t),
1574 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1575 dofstart, block_size);
1576 }
1577
1578 // Rotate a face
1579 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1580 {
1581 const auto& m = trans[0];
1582 const auto& v_size_t = std::get<0>(m);
1583 const auto& matrix = std::get<1>(m);
1584 op(std::span(v_size_t),
1585 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1586 dofstart, block_size);
1587 }
1588
1589 // Reflect a face (post rotation)
1590 if (post and cell_info >> (3 * f) & 1)
1591 {
1592 const auto& m = trans[1];
1593 const auto& v_size_t = std::get<0>(m);
1594 const auto& matrix = std::get<1>(m);
1595 op(std::span(v_size_t),
1596 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1597 dofstart, block_size);
1598 }
1599
1600 dofstart += _edofs[2][f].size();
1601 }
1602 }
1603 }
1604}
1605//-----------------------------------------------------------------------------
1606template <std::floating_point F>
1607template <typename T>
1608void FiniteElement<F>::T_apply(std::span<T> u, int n,
1609 std::uint32_t cell_info) const
1610{
1611 if (_dof_transformations_are_identity)
1612 return;
1613
1614 if (_dof_transformations_are_permutations)
1615 permute_data<T, false>(u, n, cell_info, _eperm);
1616 else
1617 {
1618 transform_data<T, false>(u, n, cell_info, _etrans,
1619 precompute::apply_matrix<F, T>);
1620 }
1621}
1622//-----------------------------------------------------------------------------
1623template <std::floating_point F>
1624template <typename T>
1625void FiniteElement<F>::Tt_apply(std::span<T> u, int n,
1626 std::uint32_t cell_info) const
1627{
1628 if (_dof_transformations_are_identity)
1629 return;
1630 else if (_dof_transformations_are_permutations)
1631 permute_data<T, true>(u, n, cell_info, _eperm_rev);
1632 else
1633 {
1634 transform_data<T, true>(u, n, cell_info, _etransT,
1635 precompute::apply_matrix<F, T>);
1636 }
1637}
1638//-----------------------------------------------------------------------------
1639template <std::floating_point F>
1640template <typename T>
1641void FiniteElement<F>::Tt_inv_apply(std::span<T> u, int n,
1642 std::uint32_t cell_info) const
1643{
1644 if (_dof_transformations_are_identity)
1645 return;
1646 else if (_dof_transformations_are_permutations)
1647 permute_data<T, false>(u, n, cell_info, _eperm);
1648 else
1649 {
1650 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1651 precompute::apply_matrix<F, T>);
1652 }
1653}
1654//-----------------------------------------------------------------------------
1655template <std::floating_point F>
1656template <typename T>
1657void FiniteElement<F>::Tinv_apply(std::span<T> u, int n,
1658 std::uint32_t cell_info) const
1659{
1660 if (_dof_transformations_are_identity)
1661 return;
1662 else if (_dof_transformations_are_permutations)
1663 permute_data<T, true>(u, n, cell_info, _eperm_rev);
1664 else
1665 {
1666 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1667 precompute::apply_matrix<F, T>);
1668 }
1669}
1670//-----------------------------------------------------------------------------
1671template <std::floating_point F>
1672template <typename T>
1673void FiniteElement<F>::Tt_apply_right(std::span<T> u, int n,
1674 std::uint32_t cell_info) const
1675{
1676 if (_dof_transformations_are_identity)
1677 return;
1678 else if (_dof_transformations_are_permutations)
1679 {
1680 assert(u.size() % n == 0);
1681 const int step = u.size() / n;
1682 for (int i = 0; i < n; ++i)
1683 {
1684 std::span<T> dblock(u.data() + i * step, step);
1685 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1686 }
1687 }
1688 else
1689 {
1690 transform_data<T, false>(u, n, cell_info, _etrans,
1691 precompute::apply_tranpose_matrix_right<F, T>);
1692 }
1693}
1694//-----------------------------------------------------------------------------
1695template <std::floating_point F>
1696template <typename T>
1697void FiniteElement<F>::Tinv_apply_right(std::span<T> u, int n,
1698 std::uint32_t cell_info) const
1699{
1700 if (_dof_transformations_are_identity)
1701 return;
1702 else if (_dof_transformations_are_permutations)
1703 {
1704 assert(u.size() % n == 0);
1705 const int step = u.size() / n;
1706 for (int i = 0; i < n; ++i)
1707 {
1708 std::span<T> dblock(u.data() + i * step, step);
1709 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1710 }
1711 }
1712 else
1713 {
1714 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1715 precompute::apply_tranpose_matrix_right<F, T>);
1716 }
1717}
1718//-----------------------------------------------------------------------------
1719template <std::floating_point F>
1720template <typename T>
1721void FiniteElement<F>::T_apply_right(std::span<T> u, int n,
1722 std::uint32_t cell_info) const
1723{
1724 if (_dof_transformations_are_identity)
1725 return;
1726 else if (_dof_transformations_are_permutations)
1727 {
1728 assert(u.size() % n == 0);
1729 const int step = u.size() / n;
1730 for (int i = 0; i < n; ++i)
1731 {
1732 std::span<T> dblock(u.data() + i * step, step);
1733 permute_data<T, true>(dblock, 1, cell_info, _eperm_rev);
1734 }
1735 }
1736 else
1737 {
1738 transform_data<T, true>(u, n, cell_info, _etransT,
1739 precompute::apply_tranpose_matrix_right<F, T>);
1740 }
1741}
1742//-----------------------------------------------------------------------------
1743template <std::floating_point F>
1744template <typename T>
1745void FiniteElement<F>::Tt_inv_apply_right(std::span<T> u, int n,
1746 std::uint32_t cell_info) const
1747{
1748 if (_dof_transformations_are_identity)
1749 return;
1750 else if (_dof_transformations_are_permutations)
1751 {
1752 assert(u.size() % n == 0);
1753 const int step = u.size() / n;
1754 for (int i = 0; i < n; ++i)
1755 {
1756 std::span<T> dblock(u.data() + i * step, step);
1757 permute_data<T, true>(dblock, 1, cell_info, _eperm_rev);
1758 }
1759 }
1760 else
1761 {
1762 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1763 precompute::apply_tranpose_matrix_right<F, T>);
1764 }
1765}
1766//-----------------------------------------------------------------------------
1767
1768} // namespace basix
A finite element.
Definition: finite-element.h:139
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.
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:1133
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:1697
std::pair< std::vector< F >, std::array< std::size_t, 3 > > base_transformations() const
Get the base transformations.
Definition: finite-element.cpp:1339
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:1608
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:1625
bool dof_transformations_are_identity() const
Indicates is the dof transformations are all the identity.
Definition: finite-element.h:557
bool operator==(const FiniteElement &e) const
Check if two elements are the same.
Definition: finite-element.cpp:1157
std::map< cell::type, std::pair< std::vector< F >, std::array< std::size_t, 3 > > > entity_transformations() const
Return the entity dof transformation matrices.
Definition: finite-element.h:768
const std::vector< int > & dof_ordering() const
Get dof layout.
Definition: finite-element.h:1192
int embedded_subdegree() const
Highest degree n such that a Lagrange (or vector Lagrange) element of degree n is a subspace of this ...
Definition: finite-element.h:505
const std::vector< std::vector< std::vector< int > > > & entity_dofs() const
Get the dofs on each topological entity: (vertices, edges, faces, cell) in that order.
Definition: finite-element.h:660
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:1721
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:1157
int interpolation_nderivs() const
The number of derivatives needed when interpolating.
Definition: finite-element.h:1189
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:1247
int embedded_superdegree() const
Lowest degree n such that the highest degree polynomial in this element is contained in a Lagrange (o...
Definition: finite-element.h:500
int dim() const
Dimension of the finite element space.
Definition: finite-element.h:518
void permute(std::span< std::int32_t > d, std::uint32_t cell_info) const
Permute indices associated with degree-of-freedoms on the reference element ordering to the globally ...
Definition: finite-element.h:795
std::size_t hash() const
Get a unique hash of this element.
Definition: finite-element.cpp:1197
void permute_inv(std::span< std::int32_t > d, std::uint32_t cell_info) const
Perform the inverse of the operation applied by permute().
Definition: finite-element.h:826
void Tt_apply_right(std::span< T > u, int n, std::uint32_t cell_info) const
Right(post)-apply the transpose of the operator applied by T_apply().
Definition: finite-element.h:1673
const std::vector< std::size_t > & value_shape() const
Element value tensor shape.
Definition: finite-element.h:512
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:1641
element::lagrange_variant lagrange_variant() const
Lagrange variant of the element.
Definition: finite-element.h:526
void Tt_inv_apply_right(std::span< T > u, int n, std::uint32_t cell_info) const
Right(post)-apply the transpose inverse of the operator applied by T_apply().
Definition: finite-element.h:1745
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:1090
int degree() const
Get the element polynomial degree.
Definition: finite-element.h:494
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:1657
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & points() const
Return the interpolation points.
Definition: finite-element.h:970
sobolev::space sobolev_space() const
Underlying Sobolev space for this element.
Definition: finite-element.h:541
const std::vector< std::vector< std::vector< int > > > & entity_closure_dofs() const
Get the dofs on the closure of each topological entity: (vertices, edges, faces, cell) in that order.
Definition: finite-element.h:674
FiniteElement(const FiniteElement &element)=default
Copy constructor.
std::array< std::size_t, 4 > tabulate_shape(std::size_t nd, std::size_t num_points) const
Array shape for tabulate basis values and derivatives at set of points.
Definition: finite-element.h:363
FiniteElement & operator=(FiniteElement &&element)=default
Move assignment operator.
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:1027
polyset::type polyset_type() const
Get the element polyset type.
Definition: finite-element.h:490
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & wcoeffs() const
Get the coefficients that define the polynomial set in terms of the orthonormal polynomials.
Definition: finite-element.h:1079
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:1448
bool dof_transformations_are_permutations() const
Indicates if the degree-of-freedom transformations are all permutations.
Definition: finite-element.h:551
std::pair< std::vector< F >, std::array< std::size_t, 3 > > push_forward(impl::mdspan_t< const F, 3 > U, impl::mdspan_t< const F, 3 > J, std::span< const F > detJ, impl::mdspan_t< const F, 3 > K) const
Map function values from the reference to a physical cell.
Definition: finite-element.cpp:1408
cell::type cell_type() const
Get the element cell type.
Definition: finite-element.h:486
bool interpolation_is_identity() const
Indicates whether or not the interpolation matrix for this element is an identity matrix.
Definition: finite-element.h:1186
bool discontinuous() const
Indicates whether this element is the discontinuous variant.
Definition: finite-element.h:547
std::vector< std::vector< FiniteElement< F > > > get_tensor_product_representation() const
Get the tensor product representation of this element.
Definition: finite-element.h:1175
maps::type map_type() const
Map type for the element.
Definition: finite-element.h:537
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & dual_matrix() const
Get the dual matrix.
Definition: finite-element.h:1038
element::dpc_variant dpc_variant() const
DPC variant of the element.
Definition: finite-element.h:533
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & coefficient_matrix() const
Get the matrix of coefficients.
Definition: finite-element.h:1142
FiniteElement & operator=(const FiniteElement &element)=default
Assignment operator.
element::family family() const
The finite element family.
Definition: finite-element.h:522
~FiniteElement()=default
Destructor.
std::function< void(O &, const P &, const Q &, F, const R &)> map_fn() const
Return a function that performs the appropriate push-forward/pull-back for the element type.
Definition: finite-element.h:623
type
Cell type.
Definition: cell.h:21
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:520
lagrange_variant
Variants of a Lagrange space that can be created.
Definition: element-families.h:12
impl::mdspan_t< T, d > mdspan_t
Typedef for mdspan.
Definition: finite-element.h:104
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:60
void contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Contravariant Piola map.
Definition: maps.h:80
void double_contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Double contravariant Piola map.
Definition: maps.h:134
void double_covariant_piola(O &&r, const P &U, const Q &J, double, const R &K)
Double covariant Piola map.
Definition: maps.h:102
type
Map type.
Definition: maps.h:38
type
Cell type.
Definition: polyset.h:136
void apply_permutation_mapped(std::span< const std::size_t > perm, std::span< E > data, std::span< const int > emap, std::size_t n=1)
Permutation of mapped data.
Definition: precompute.h:154
space
Sobolev space type.
Definition: sobolev-spaces.h:13
Basix: FEniCS runtime basis evaluation library.
Definition: cell.h:17
FiniteElement< T > create_element(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous, std::vector< int > dof_ordering={})
Definition: finite-element.cpp:193
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:399
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:349
std::string version()
Definition: finite-element.cpp:1486
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:609
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:330