Basix
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{
29namespace impl
30{
31template <typename T, std::size_t d>
32using mdspan_t = md::mdspan<T, md::dextents<std::size_t, d>>;
33template <typename T, std::size_t d>
34using mdarray_t
35 = md::MDSPAN_IMPL_PROPOSED_NAMESPACE::mdarray<T,
36 md::dextents<std::size_t, d>>;
37
40template <typename T>
41std::array<std::vector<mdspan_t<const T, 2>>, 4>
42to_mdspan(std::array<std::vector<mdarray_t<T, 2>>, 4>& x)
43{
44 std::array<std::vector<mdspan_t<const T, 2>>, 4> x1;
45 for (std::size_t i = 0; i < x.size(); ++i)
46 for (std::size_t j = 0; j < x[i].size(); ++j)
47 x1[i].emplace_back(x[i][j].data(), x[i][j].extents());
48
49 return x1;
50}
51
54template <typename T>
55std::array<std::vector<mdspan_t<const T, 4>>, 4>
56to_mdspan(std::array<std::vector<mdarray_t<T, 4>>, 4>& M)
57{
58 std::array<std::vector<mdspan_t<const T, 4>>, 4> M1;
59 for (std::size_t i = 0; i < M.size(); ++i)
60 for (std::size_t j = 0; j < M[i].size(); ++j)
61 M1[i].emplace_back(M[i][j].data(), M[i][j].extents());
62
63 return M1;
64}
65
68template <typename T>
69std::array<std::vector<mdspan_t<const T, 2>>, 4>
70to_mdspan(const std::array<std::vector<std::vector<T>>, 4>& x,
71 const std::array<std::vector<std::array<std::size_t, 2>>, 4>& shape)
72{
73 std::array<std::vector<mdspan_t<const T, 2>>, 4> x1;
74 for (std::size_t i = 0; i < x.size(); ++i)
75 for (std::size_t j = 0; j < x[i].size(); ++j)
76 x1[i].push_back(mdspan_t<const T, 2>(x[i][j].data(), shape[i][j]));
77
78 return x1;
79}
80
83template <typename T>
84std::array<std::vector<mdspan_t<const T, 4>>, 4>
85to_mdspan(const std::array<std::vector<std::vector<T>>, 4>& M,
86 const std::array<std::vector<std::array<std::size_t, 4>>, 4>& shape)
87{
88 std::array<std::vector<mdspan_t<const T, 4>>, 4> M1;
89 for (std::size_t i = 0; i < M.size(); ++i)
90 for (std::size_t j = 0; j < M[i].size(); ++j)
91 M1[i].push_back(mdspan_t<const T, 4>(M[i][j].data(), shape[i][j]));
92
93 return M1;
94}
95
96} // namespace impl
97
98namespace element
99{
101template <typename T, std::size_t d>
102using mdspan_t = impl::mdspan_t<T, d>;
103
119template <std::floating_point T>
120std::tuple<std::array<std::vector<std::vector<T>>, 4>,
121 std::array<std::vector<std::array<std::size_t, 2>>, 4>,
122 std::array<std::vector<std::vector<T>>, 4>,
123 std::array<std::vector<std::array<std::size_t, 4>>, 4>>
124make_discontinuous(const std::array<std::vector<mdspan_t<const T, 2>>, 4>& x,
125 const std::array<std::vector<mdspan_t<const T, 4>>, 4>& M,
126 std::size_t tdim, std::size_t value_size);
127
128} // namespace element
129
135template <std::floating_point F>
137{
138 template <typename T, std::size_t d>
139 using mdspan_t = md::mdspan<T, md::dextents<std::size_t, d>>;
140
141public:
143 using scalar_type = F;
144
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
332
335
337 ~FiniteElement() = default;
338
341
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_inv);
838 }
839
855 void permute_subentity_closure(std::span<std::int32_t> d,
856 std::uint32_t cell_info,
857 cell::type entity_type, int entity_index) const
858 {
859 const int entity_dim = cell::topological_dimension(entity_type);
860
861 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
862
863 std::uint32_t entity_info = 0;
864 switch (entity_dim)
865 {
866 case 0:
867 return;
868 case 1:
869 entity_info = cell_info >> (face_start + entity_index) & 1;
870 break;
871 case 2:
872 entity_info = cell_info >> (3 * entity_index) & 7;
873 break;
874 default:
875 throw std::runtime_error("Unsupported cell dimension");
876 }
877 permute_subentity_closure(d, entity_info, entity_type);
878 }
879
892 void permute_subentity_closure_inv(std::span<std::int32_t> d,
893 std::uint32_t cell_info,
894 cell::type entity_type,
895 int entity_index) const
896 {
897 const int entity_dim = cell::topological_dimension(entity_type);
898
899 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
900
901 std::uint32_t entity_info;
902 switch (entity_dim)
903 {
904 case 0:
905 return;
906 case 1:
907 entity_info = cell_info >> (face_start + entity_index) & 1;
908 break;
909 case 2:
910 entity_info = cell_info >> (3 * entity_index) & 7;
911 break;
912 default:
913 throw std::runtime_error("Unsupported cell dimension");
914 }
915 permute_subentity_closure_inv(d, entity_info, entity_type);
916 }
917
932
933 void permute_subentity_closure(std::span<std::int32_t> d,
934 std::uint32_t entity_info,
935 cell::type entity_type) const
936 {
937 if (!_dof_transformations_are_permutations)
938 {
939 throw std::runtime_error(
940 "The DOF transformations for this element are not permutations");
941 }
942
943 const int entity_dim = cell::topological_dimension(entity_type);
944
945 if (entity_dim == 0)
946 return;
947
948 auto& perm = _subentity_closure_perm.at(entity_type);
949 if (entity_dim == 1)
950 {
951 if (entity_info & 1)
952 {
954 }
955 }
956 else if (entity_dim == 2)
957 {
958 // Rotate a face
959 for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
960 {
962 }
963
964 // Reflect a face (post rotate)
965 if (entity_info & 1)
966 {
968 }
969 }
970 else
971 {
972 throw std::runtime_error(
973 "Invalid dimension for permute_subentity_closure");
974 }
975 }
976
988 void permute_subentity_closure_inv(std::span<std::int32_t> d,
989 std::uint32_t entity_info,
990 cell::type entity_type) const
991 {
992 if (!_dof_transformations_are_permutations)
993 {
994 throw std::runtime_error(
995 "The DOF transformations for this element are not permutations");
996 }
997
998 const int entity_dim = cell::topological_dimension(entity_type);
999
1000 if (entity_dim == 0)
1001 return;
1002
1003 auto& perm = _subentity_closure_perm_inv.at(entity_type);
1004 if (entity_dim == 1)
1005 {
1006 if (entity_info & 1)
1007 {
1009 }
1010 }
1011 else if (entity_dim == 2)
1012 {
1013 // Reflect a face (pre rotate)
1014 if (entity_info & 1)
1015 {
1017 }
1018
1019 // Rotate a face
1020 for (std::uint32_t r = 0; r < (entity_info >> 1 & 3); ++r)
1021 {
1023 }
1024 }
1025 else
1026 {
1027 throw std::runtime_error(
1028 "Invalid dimension for permute_subentity_closure");
1029 }
1030 }
1031
1062 template <typename T>
1063 void T_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
1064
1075 template <typename T>
1076 void Tt_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
1077
1089 template <typename T>
1090 void Tt_inv_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
1091
1102 template <typename T>
1103 void Tinv_apply(std::span<T> u, int n, std::uint32_t cell_info) const;
1104
1115 template <typename T>
1116 void Tt_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
1117
1127 template <typename T>
1128 void T_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
1129
1140 template <typename T>
1141 void Tinv_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
1142
1153 template <typename T>
1154 void Tt_inv_apply_right(std::span<T> u, int n, std::uint32_t cell_info) const;
1155
1162 const std::pair<std::vector<F>, std::array<std::size_t, 2>>& points() const
1163 {
1164 return _points;
1165 }
1166
1218 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1220 {
1221 return _matM;
1222 }
1223
1229 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1231 {
1232 return _dual_matrix;
1233 }
1234
1271 const std::pair<std::vector<F>, std::array<std::size_t, 2>>& wcoeffs() const
1272 {
1273 return _wcoeffs;
1274 }
1275
1280 const std::array<
1281 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>, 4>&
1282 x() const
1283 {
1284 return _x;
1285 }
1286
1323 const std::array<
1324 std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>, 4>&
1325 M() const
1326 {
1327 return _M;
1328 }
1329
1333 const std::pair<std::vector<F>, std::array<std::size_t, 2>>&
1335 {
1336 return _coeffs;
1337 }
1338
1350 {
1351 return !_tensor_factors.empty();
1352 }
1353
1366 std::vector<std::vector<FiniteElement<F>>>
1368 {
1370 throw std::runtime_error("Element has no tensor product representation.");
1371 return _tensor_factors;
1372 }
1373
1378 bool interpolation_is_identity() const { return _interpolation_is_identity; }
1379
1381 int interpolation_nderivs() const { return _interpolation_nderivs; }
1382
1384 const std::vector<int>& dof_ordering() const { return _dof_ordering; }
1385
1386private:
1393 template <typename T, bool post>
1394 void permute_data(
1395 std::span<T> data, int block_size, std::uint32_t cell_info,
1396 const std::map<cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1397 const;
1398
1399 using array2_t = std::pair<std::vector<F>, std::array<std::size_t, 2>>;
1400 using array3_t = std::pair<std::vector<F>, std::array<std::size_t, 3>>;
1401 using trans_data_t
1402 = std::vector<std::pair<std::vector<std::size_t>, array2_t>>;
1403
1410 template <typename T, bool post, typename OP>
1411 void
1412 transform_data(std::span<T> data, int block_size, std::uint32_t cell_info,
1413 const std::map<cell::type, trans_data_t>& etrans, OP op) const;
1414
1415 // Cell type
1416 cell::type _cell_type;
1417
1418 // Polyset type
1419 polyset::type _poly_type;
1420
1421 // Topological dimension of the cell
1422 std::size_t _cell_tdim;
1423
1424 // Topological dimension of the cell
1425 std::vector<std::vector<cell::type>> _cell_subentity_types;
1426
1427 // Finite element family
1428 element::family _family;
1429
1430 // Lagrange variant
1431 element::lagrange_variant _lagrange_variant;
1432
1433 // DPC variant
1434 element::dpc_variant _dpc_variant;
1435
1436 // Degree that was input when creating the element
1437 int _degree;
1438
1439 // Degree
1440 int _interpolation_nderivs;
1441
1442 // Highest degree polynomial in element's polyset
1443 int _embedded_superdegree;
1444
1445 // Highest degree space that is a subspace of element's polyset
1446 int _embedded_subdegree;
1447
1448 // Value shape
1449 std::vector<std::size_t> _value_shape;
1450
1452 maps::type _map_type;
1453
1455 sobolev::space _sobolev_space;
1456
1457 // Shape function coefficient of expansion sets on cell. If shape
1458 // function is given by @f$\psi_i = \sum_{k} \phi_{k}
1459 // \alpha^{i}_{k}@f$, then _coeffs(i, j) = @f$\alpha^i_k@f$. ie
1460 // _coeffs.row(i) are the expansion coefficients for shape function i
1461 // (@f$\psi_{i}@f$).
1462 std::pair<std::vector<F>, std::array<std::size_t, 2>> _coeffs;
1463
1464 // Dofs associated with each cell (sub-)entity
1465 std::vector<std::vector<std::vector<int>>> _edofs;
1466
1467 // Dofs associated with the closdure of each cell (sub-)entity
1468 std::vector<std::vector<std::vector<int>>> _e_closure_dofs;
1469
1470 // Entity transformations
1471 std::map<cell::type, array3_t> _entity_transformations;
1472
1473 // Set of points used for point evaluation
1474 // Experimental - currently used for an implementation of
1475 // "tabulate_dof_coordinates" Most useful for Lagrange. This may change or go
1476 // away. For non-Lagrange elements, these points will be used in combination
1477 // with _interpolation_matrix to perform interpolation
1478 std::pair<std::vector<F>, std::array<std::size_t, 2>> _points;
1479
1480 // Interpolation points on the cell. The shape is (entity_dim, num
1481 // entities of given dimension, num_points, tdim)
1482 std::array<std::vector<std::pair<std::vector<F>, std::array<std::size_t, 2>>>,
1483 4>
1484 _x;
1485
1487 std::pair<std::vector<F>, std::array<std::size_t, 2>> _matM;
1488
1489 // Indicates whether or not the DOF transformations are all
1490 // permutations
1491 bool _dof_transformations_are_permutations;
1492
1493 // Indicates whether or not the DOF transformations are all identity
1494 bool _dof_transformations_are_identity;
1495
1496 // The entity permutations (factorised). This will only be set if
1497 // _dof_transformations_are_permutations is True and
1498 // _dof_transformations_are_identity is False
1499 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm;
1500
1501 // The reverse entity permutations (factorised). This will only be set
1502 // if _dof_transformations_are_permutations is True and
1503 // _dof_transformations_are_identity is False
1504 std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm_inv;
1505
1506 // The entity transformations in precomputed form
1507 std::map<cell::type, trans_data_t> _etrans;
1508
1509 // The transposed entity transformations in precomputed form
1510 std::map<cell::type, trans_data_t> _etransT;
1511
1512 // The inverse entity transformations in precomputed form
1513 std::map<cell::type, trans_data_t> _etrans_inv;
1514
1515 // The inverse transpose entity transformations in precomputed form
1516 std::map<cell::type, trans_data_t> _etrans_invT;
1517
1518 // The subentity closure permutations (factorised). This will only be set if
1519 // _dof_transformations_are_permutations is True
1520 std::map<cell::type, std::vector<std::vector<std::size_t>>>
1521 _subentity_closure_perm;
1522
1523 // The inverse subentity closure permutations (factorised). This will only be
1524 // set if _dof_transformations_are_permutations is True
1525 std::map<cell::type, std::vector<std::vector<std::size_t>>>
1526 _subentity_closure_perm_inv;
1527
1528 // Indicates whether or not this is the discontinuous version of the
1529 // element
1530 bool _discontinuous;
1531
1532 // The dual matrix
1533 std::pair<std::vector<F>, std::array<std::size_t, 2>> _dual_matrix;
1534
1535 // Dof reordering for different element dof layout compatibility.
1536 // The reference basix layout is ordered by entity, i.e. dofs on
1537 // vertices, followed by edges, faces, then internal dofs.
1538 // _dof_ordering stores the map to the new order required, e.g.
1539 // for a P2 triangle, _dof_ordering=[0 3 5 1 2 4] will place
1540 // dofs 0, 3, 5 on the vertices and 1, 2, 4, on the edges.
1541 std::vector<int> _dof_ordering;
1542
1543 // Tensor product representation
1544 // Entries of tuple are (list of elements on an interval, permutation
1545 // of DOF numbers)
1546 // @todo: For vector-valued elements, a tensor product type and a
1547 // scaling factor may additionally be needed.
1548 std::vector<std::vector<FiniteElement>> _tensor_factors;
1549
1550 // Is the interpolation matrix an identity?
1551 bool _interpolation_is_identity;
1552
1553 // The coefficients that define the polynomial set in terms of the
1554 // orthonormal polynomials
1555 std::pair<std::vector<F>, std::array<std::size_t, 2>> _wcoeffs;
1556
1557 // Interpolation matrices for each entity
1558 using array4_t
1559 = std::vector<std::pair<std::vector<F>, std::array<std::size_t, 4>>>;
1560 std::array<array4_t, 4> _M;
1561};
1562
1588template <std::floating_point T>
1589FiniteElement<T> create_custom_element(
1590 cell::type cell_type, const std::vector<std::size_t>& value_shape,
1591 impl::mdspan_t<const T, 2> wcoeffs,
1592 const std::array<std::vector<impl::mdspan_t<const T, 2>>, 4>& x,
1593 const std::array<std::vector<impl::mdspan_t<const T, 4>>, 4>& M,
1594 int interpolation_nderivs, maps::type map_type,
1595 sobolev::space sobolev_space, bool discontinuous, int embedded_subdegree,
1596 int embedded_superdegree, polyset::type poly_type);
1597
1609template <std::floating_point T>
1610FiniteElement<T> create_element(element::family family, cell::type cell,
1611 int degree, element::lagrange_variant lvariant,
1612 element::dpc_variant dvariant,
1613 bool discontinuous,
1614 std::vector<int> dof_ordering = {});
1615
1627std::vector<int> tp_dof_ordering(element::family family, cell::type cell,
1628 int degree, element::lagrange_variant lvariant,
1629 element::dpc_variant dvariant,
1630 bool discontinuous);
1631
1642std::vector<int> lex_dof_ordering(element::family family, cell::type cell,
1643 int degree, element::lagrange_variant lvariant,
1644 element::dpc_variant dvariant,
1645 bool discontinuous);
1646
1659template <std::floating_point T>
1660std::vector<std::vector<FiniteElement<T>>>
1661tp_factors(element::family family, cell::type cell, int degree,
1663 bool discontinuous, const std::vector<int>& dof_ordering);
1664
1676template <std::floating_point T>
1680 element::dpc_variant dvariant, bool discontinuous);
1681
1684std::string version();
1685
1686//-----------------------------------------------------------------------------
1687template <std::floating_point F>
1688template <typename T, bool post>
1689void FiniteElement<F>::permute_data(
1690 std::span<T> data, int block_size, std::uint32_t cell_info,
1691 const std::map<cell::type, std::vector<std::vector<std::size_t>>>& eperm)
1692 const
1693{
1694 if (_cell_tdim >= 2)
1695 {
1696 // This assumes 3 bits are used per face. This will need updating if 3D
1697 // cells with faces with more than 4 sides are implemented
1698 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1699
1700 // Permute DOFs on edges
1701 {
1702 auto& trans = eperm.at(cell::type::interval)[0];
1703 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1704 {
1705 // Reverse an edge
1706 if (cell_info >> (face_start + e) & 1)
1707 {
1708 precompute::apply_permutation_mapped(trans, data, _edofs[1][e],
1709 block_size);
1710 }
1711 }
1712 }
1713
1714 if (_cell_tdim == 3)
1715 {
1716 // Permute DOFs on faces
1717 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1718 {
1719 auto& trans = eperm.at(_cell_subentity_types[2][f]);
1720
1721 // Reflect a face (pre rotate)
1722 if (!post and cell_info >> (3 * f) & 1)
1723 {
1724 precompute::apply_permutation_mapped(trans[1], data, _edofs[2][f],
1725 block_size);
1726 }
1727
1728 // Rotate a face
1729 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1730 {
1731 precompute::apply_permutation_mapped(trans[0], data, _edofs[2][f],
1732 block_size);
1733 }
1734
1735 // Reflect a face (post rotate)
1736 if (post and cell_info >> (3 * f) & 1)
1737 {
1738 precompute::apply_permutation_mapped(trans[1], data, _edofs[2][f],
1739 block_size);
1740 }
1741 }
1742 }
1743 }
1744}
1745//-----------------------------------------------------------------------------
1746template <std::floating_point F>
1747template <typename T, bool post, typename OP>
1748void FiniteElement<F>::transform_data(
1749 std::span<T> data, int block_size, std::uint32_t cell_info,
1750 const std::map<cell::type, trans_data_t>& etrans, OP op) const
1751{
1752 if (_cell_tdim >= 2)
1753 {
1754 // This assumes 3 bits are used per face. This will need updating if
1755 // 3D cells with faces with more than 4 sides are implemented
1756 int face_start = _cell_tdim == 3 ? 3 * _edofs[2].size() : 0;
1757 int dofstart = 0;
1758 for (auto& edofs0 : _edofs[0])
1759 dofstart += edofs0.size();
1760
1761 // Transform DOFs on edges
1762 {
1763 auto& [v_size_t, matrix] = etrans.at(cell::type::interval)[0];
1764 for (std::size_t e = 0; e < _edofs[1].size(); ++e)
1765 {
1766 // Reverse an edge
1767 if (cell_info >> (face_start + e) & 1)
1768 {
1769 op(std::span(v_size_t),
1770 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1771 dofstart, block_size);
1772 }
1773 dofstart += _edofs[1][e].size();
1774 }
1775 }
1776
1777 if (_cell_tdim == 3)
1778 {
1779 // Permute DOFs on faces
1780 for (std::size_t f = 0; f < _edofs[2].size(); ++f)
1781 {
1782 auto& trans = etrans.at(_cell_subentity_types[2][f]);
1783
1784 // Reflect a face (pre rotation)
1785 if (!post and cell_info >> (3 * f) & 1)
1786 {
1787 const auto& m = trans[1];
1788 const auto& v_size_t = std::get<0>(m);
1789 const auto& matrix = std::get<1>(m);
1790 op(std::span(v_size_t),
1791 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1792 dofstart, block_size);
1793 }
1794
1795 // Rotate a face
1796 for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1797 {
1798 const auto& m = trans[0];
1799 const auto& v_size_t = std::get<0>(m);
1800 const auto& matrix = std::get<1>(m);
1801 op(std::span(v_size_t),
1802 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1803 dofstart, block_size);
1804 }
1805
1806 // Reflect a face (post rotation)
1807 if (post and cell_info >> (3 * f) & 1)
1808 {
1809 const auto& m = trans[1];
1810 const auto& v_size_t = std::get<0>(m);
1811 const auto& matrix = std::get<1>(m);
1812 op(std::span(v_size_t),
1813 mdspan_t<const F, 2>(matrix.first.data(), matrix.second), data,
1814 dofstart, block_size);
1815 }
1816
1817 dofstart += _edofs[2][f].size();
1818 }
1819 }
1820 }
1821}
1822//-----------------------------------------------------------------------------
1823template <std::floating_point F>
1824template <typename T>
1825void FiniteElement<F>::T_apply(std::span<T> u, int n,
1826 std::uint32_t cell_info) const
1827{
1828 if (_dof_transformations_are_identity)
1829 return;
1830
1831 if (_dof_transformations_are_permutations)
1832 permute_data<T, false>(u, n, cell_info, _eperm);
1833 else
1834 {
1835 transform_data<T, false>(u, n, cell_info, _etrans,
1837 }
1838}
1839//-----------------------------------------------------------------------------
1840template <std::floating_point F>
1841template <typename T>
1842void FiniteElement<F>::Tt_apply(std::span<T> u, int n,
1843 std::uint32_t cell_info) const
1844{
1845 if (_dof_transformations_are_identity)
1846 return;
1847 else if (_dof_transformations_are_permutations)
1848 permute_data<T, true>(u, n, cell_info, _eperm_inv);
1849 else
1850 {
1851 transform_data<T, true>(u, n, cell_info, _etransT,
1853 }
1854}
1855//-----------------------------------------------------------------------------
1856template <std::floating_point F>
1857template <typename T>
1858void FiniteElement<F>::Tt_inv_apply(std::span<T> u, int n,
1859 std::uint32_t cell_info) const
1860{
1861 if (_dof_transformations_are_identity)
1862 return;
1863 else if (_dof_transformations_are_permutations)
1864 permute_data<T, false>(u, n, cell_info, _eperm);
1865 else
1866 {
1867 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1869 }
1870}
1871//-----------------------------------------------------------------------------
1872template <std::floating_point F>
1873template <typename T>
1874void FiniteElement<F>::Tinv_apply(std::span<T> u, int n,
1875 std::uint32_t cell_info) const
1876{
1877 if (_dof_transformations_are_identity)
1878 return;
1879 else if (_dof_transformations_are_permutations)
1880 permute_data<T, true>(u, n, cell_info, _eperm_inv);
1881 else
1882 {
1883 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1885 }
1886}
1887//-----------------------------------------------------------------------------
1888template <std::floating_point F>
1889template <typename T>
1890void FiniteElement<F>::Tt_apply_right(std::span<T> u, int n,
1891 std::uint32_t cell_info) const
1892{
1893 if (_dof_transformations_are_identity)
1894 return;
1895 else if (_dof_transformations_are_permutations)
1896 {
1897 assert(u.size() % n == 0);
1898 const int step = u.size() / n;
1899 for (int i = 0; i < n; ++i)
1900 {
1901 std::span<T> dblock(u.data() + i * step, step);
1902 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1903 }
1904 }
1905 else
1906 {
1907 transform_data<T, false>(u, n, cell_info, _etrans,
1909 }
1910}
1911//-----------------------------------------------------------------------------
1912template <std::floating_point F>
1913template <typename T>
1914void FiniteElement<F>::Tinv_apply_right(std::span<T> u, int n,
1915 std::uint32_t cell_info) const
1916{
1917 if (_dof_transformations_are_identity)
1918 return;
1919 else if (_dof_transformations_are_permutations)
1920 {
1921 assert(u.size() % n == 0);
1922 const int step = u.size() / n;
1923 for (int i = 0; i < n; ++i)
1924 {
1925 std::span<T> dblock(u.data() + i * step, step);
1926 permute_data<T, false>(dblock, 1, cell_info, _eperm);
1927 }
1928 }
1929 else
1930 {
1931 transform_data<T, false>(u, n, cell_info, _etrans_invT,
1933 }
1934}
1935//-----------------------------------------------------------------------------
1936template <std::floating_point F>
1937template <typename T>
1938void FiniteElement<F>::T_apply_right(std::span<T> u, int n,
1939 std::uint32_t cell_info) const
1940{
1941 if (_dof_transformations_are_identity)
1942 return;
1943 else if (_dof_transformations_are_permutations)
1944 {
1945 assert(u.size() % n == 0);
1946 const int step = u.size() / n;
1947 for (int i = 0; i < n; ++i)
1948 {
1949 std::span<T> dblock(u.data() + i * step, step);
1950 permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1951 }
1952 }
1953 else
1954 {
1955 transform_data<T, true>(u, n, cell_info, _etransT,
1957 }
1958}
1959//-----------------------------------------------------------------------------
1960template <std::floating_point F>
1961template <typename T>
1962void FiniteElement<F>::Tt_inv_apply_right(std::span<T> u, int n,
1963 std::uint32_t cell_info) const
1964{
1965 if (_dof_transformations_are_identity)
1966 return;
1967 else if (_dof_transformations_are_permutations)
1968 {
1969 assert(u.size() % n == 0);
1970 const int step = u.size() / n;
1971 for (int i = 0; i < n; ++i)
1972 {
1973 std::span<T> dblock(u.data() + i * step, step);
1974 permute_data<T, true>(dblock, 1, cell_info, _eperm_inv);
1975 }
1976 }
1977 else
1978 {
1979 transform_data<T, true>(u, n, cell_info, _etrans_inv,
1981 }
1982}
1983//-----------------------------------------------------------------------------
1984
1985} // namespace basix
A finite element.
Definition finite-element.h:137
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:1325
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:1914
std::pair< std::vector< F >, std::array< std::size_t, 3 > > base_transformations() const
Get the base transformations.
Definition finite-element.cpp:1888
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:1825
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:1842
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:892
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:1710
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:1384
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:988
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:1938
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:1349
int interpolation_nderivs() const
The number of derivatives needed when interpolating.
Definition finite-element.h:1381
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:1796
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:1749
F scalar_type
Scalar type.
Definition finite-element.h:143
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:1890
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:1858
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:855
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:1962
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:1282
int degree() const
Get the element polynomial degree.
Definition finite-element.h:494
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:933
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:1874
const std::pair< std::vector< F >, std::array< std::size_t, 2 > > & points() const
Return the interpolation points.
Definition finite-element.h:1162
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:1219
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:1271
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:1993
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:1957
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:1378
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:1367
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:1230
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:1334
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
Information about reference cells.
Definition cell.h:17
type
Cell type.
Definition cell.h:21
int topological_dimension(cell::type celltype)
Definition cell.cpp:294
Interfaces for creating finite elements.
Definition e-brezzi-douglas-marini.h:13
impl::mdspan_t< T, d > mdspan_t
Typedef for mdspan.
Definition finite-element.h:102
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:819
lagrange_variant
Variants of a Lagrange space that can be created.
Definition element-families.h:12
dpc_variant
Definition element-families.h:32
family
Available element families.
Definition element-families.h:45
void covariant_piola(O &&r, const P &U, const Q &, double, const R &K)
Covariant Piola map.
Definition maps.h:62
void contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Contravariant Piola map.
Definition maps.h:82
void double_contravariant_piola(O &&r, const P &U, const Q &J, double detJ, const R &)
Double contravariant Piola map.
Definition maps.h:132
void double_covariant_piola(O &&r, const P &U, const Q &, double, const R &K)
Double covariant Piola map.
Definition maps.h:104
type
Map type.
Definition maps.h:40
type
Cell type.
Definition polyset.h:137
void apply_permutation(std::span< const std::size_t > perm, std::span< E > data, std::size_t offset=0, std::size_t n=1)
Apply a (precomputed) permutation .
Definition precompute.h:137
void apply_permutation_mapped(std::span< const std::size_t > perm, std::span< E > data, std::span< const int > emap, std::size_t n=1)
Permutation of mapped data.
Definition precompute.h:147
void apply_tranpose_matrix_right(std::span< const std::size_t > v_size_t, md::mdspan< const T, md::dextents< std::size_t, 2 > > M, std::span< E > data, std::size_t offset=0, std::size_t n=1)
Apply a (precomputed) matrix to some transposed data.
Definition precompute.h:276
void apply_matrix(std::span< const std::size_t > v_size_t, md::mdspan< const T, md::dextents< std::size_t, 2 > > M, std::span< E > data, std::size_t offset=0, std::size_t n=1)
Apply a (precomputed) matrix.
Definition precompute.h:236
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:190
std::vector< std::vector< FiniteElement< T > > > tp_factors(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous, const std::vector< int > &dof_ordering)
Definition finite-element.cpp:346
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:396
std::string version()
Definition finite-element.cpp:2027
std::vector< int > lex_dof_ordering(element::family family, cell::type cell, int degree, element::lagrange_variant lvariant, element::dpc_variant dvariant, bool discontinuous)
Definition finite-element.cpp:512
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:904
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:327