Home Installation C++ docs Python docs
finite-element.h
1 // Copyright (c) 2020 Chris Richardson
2 // FEniCS Project
3 // SPDX-License-Identifier: MIT
4 
5 // FIXME: just include everything for now
6 // Need to define public API
7 
8 #pragma once
9 
10 #include "cell.h"
11 #include "element-families.h"
12 #include "maps.h"
13 #include "precompute.h"
14 #include <array>
15 #include <numeric>
16 #include <set>
17 #include <string>
18 #include <vector>
19 #include <xtensor/xadapt.hpp>
20 #include <xtensor/xcomplex.hpp>
21 #include <xtensor/xtensor.hpp>
22 #include <xtensor/xview.hpp>
23 #include <xtl/xspan.hpp>
24 
26 namespace basix
27 {
28 
184 xt::xtensor<double, 3> compute_expansion_coefficients(
185  cell::type cell_type, const xt::xtensor<double, 2>& B,
186  const std::vector<std::vector<xt::xtensor<double, 3>>>& M,
187  const std::vector<std::vector<xt::xtensor<double, 2>>>& x, int degree,
188  double kappa_tol = 0.0);
189 
194 {
195 
196 public:
211  FiniteElement(element::family family, cell::type cell_type, int degree,
212  const std::vector<std::size_t>& value_shape,
213  const xt::xtensor<double, 3>& coeffs,
214  const std::map<cell::type, xt::xtensor<double, 3>>&
216  const std::array<std::vector<xt::xtensor<double, 2>>, 4>& x,
217  const std::array<std::vector<xt::xtensor<double, 3>>, 4>& M,
218  maps::type map_type = maps::type::identity);
219 
221  FiniteElement(const FiniteElement& element) = default;
222 
224  FiniteElement(FiniteElement&& element) = default;
225 
227  ~FiniteElement() = default;
228 
230  FiniteElement& operator=(const FiniteElement& element) = default;
231 
233  FiniteElement& operator=(FiniteElement&& element) = default;
234 
252  xt::xtensor<double, 4> tabulate(int nd, const xt::xarray<double>& x) const;
253 
258  void tabulate(int nd, const xt::xarray<double>& x,
259  xt::xtensor<double, 4>& basis_data) const;
260 
263  cell::type cell_type() const;
264 
267  int degree() const;
268 
272  int value_size() const;
273 
276  const std::vector<int>& value_shape() const;
277 
281  int dim() const;
282 
285  element::family family() const;
286 
289  maps::type mapping_type() const;
290 
294 
298 
313  xt::xtensor<double, 3>
314  map_push_forward(const xt::xtensor<double, 3>& U,
315  const xt::xtensor<double, 3>& J,
316  const xtl::span<const double>& detJ,
317  const xt::xtensor<double, 3>& K) const;
318 
334  template <typename O, typename P, typename Q, typename S, typename T>
335  void map_push_forward_m(const O& U, const P& J, const Q& detJ, const S& K,
336  T&& u) const
337  {
338  // FIXME: Should U.shape(2) be replaced by the physical value size?
339  // Can it differ?
340 
341  // Loop over points that share J
342  for (std::size_t i = 0; i < U.shape(0); ++i)
343  {
344  auto _J = xt::view(J, i, xt::all(), xt::all());
345  auto _K = xt::view(K, i, xt::all(), xt::all());
346  auto _U = xt::view(U, i, xt::all(), xt::all());
347  auto _u = xt::view(u, i, xt::all(), xt::all());
348  maps::apply_map(_u, _U, _J, detJ[i], _K, map_type);
349  }
350  }
351 
358  xt::xtensor<double, 3> map_pull_back(const xt::xtensor<double, 3>& u,
359  const xt::xtensor<double, 3>& J,
360  const xtl::span<const double>& detJ,
361  const xt::xtensor<double, 3>& K) const;
362 
378  template <typename O, typename P, typename Q, typename S, typename T>
379  void map_pull_back_m(const O& u, const P& J, const Q& detJ, const S& K,
380  T&& U) const
381  {
382  // Loop over points that share K and K
383  for (std::size_t i = 0; i < u.shape(0); ++i)
384  {
385  auto _J = xt::view(J, i, xt::all(), xt::all());
386  auto _K = xt::view(K, i, xt::all(), xt::all());
387  auto _u = xt::view(u, i, xt::all(), xt::all());
388  auto _U = xt::view(U, i, xt::all(), xt::all());
389  maps::apply_map(_U, _u, _K, 1.0 / detJ[i], _J, map_type);
390  }
391  }
392 
405  const std::vector<std::vector<int>>& num_entity_dofs() const;
406 
412  const std::vector<std::vector<int>>& num_entity_closure_dofs() const;
413 
420  const std::vector<std::vector<std::set<int>>>& entity_dofs() const;
421 
428  const std::vector<std::vector<std::set<int>>>& entity_closure_dofs() const;
429 
508  xt::xtensor<double, 3> base_transformations() const;
509 
511  std::map<cell::type, xt::xtensor<double, 3>> entity_transformations() const;
512 
516  void permute_dofs(const xtl::span<std::int32_t>& dofs,
517  std::uint32_t cell_info) const;
518 
522  void unpermute_dofs(const xtl::span<std::int32_t>& dofs,
523  std::uint32_t cell_info) const;
524 
529  template <typename T>
530  void apply_dof_transformation(const xtl::span<T>& data, int block_size,
531  std::uint32_t cell_info) const;
532 
537  template <typename T>
538  void apply_transpose_dof_transformation(const xtl::span<T>& data,
539  int block_size,
540  std::uint32_t cell_info) const;
541 
546  template <typename T>
548  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const;
549 
554  template <typename T>
555  void apply_inverse_dof_transformation(const xtl::span<T>& data,
556  int block_size,
557  std::uint32_t cell_info) const;
558 
563  template <typename T>
564  void apply_dof_transformation_to_transpose(const xtl::span<T>& data,
565  int block_size,
566  std::uint32_t cell_info) const;
567 
572  template <typename T>
574  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const;
575 
580  template <typename T>
582  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const;
583 
588  template <typename T>
590  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const;
591 
596  const xt::xtensor<double, 2>& points() const;
597 
599  int num_points() const;
600 
607  const xt::xtensor<double, 2>& interpolation_matrix() const;
608 
611 
612 private:
613  // Cell type
614  cell::type _cell_type;
615 
616  // Topological dimension of the cell
617  std::size_t _cell_tdim;
618 
619  // Topological dimension of the cell
620  std::vector<std::vector<cell::type>> _cell_subentity_types;
621 
622  // Finite element family
623  element::family _family;
624 
625  // Degree
626  int _degree;
627 
628  // Value shape
629  std::vector<int> _value_shape;
630 
632  maps::type _map_type;
633 
634  // Shape function coefficient of expansion sets on cell. If shape
635  // function is given by @f$\psi_i = \sum_{k} \phi_{k}
636  // \alpha^{i}_{k}@f$, then _coeffs(i, j) = @f$\alpha^i_k@f$. i.e.,
637  // _coeffs.row(i) are the expansion coefficients for shape function i
638  // (@f$\psi_{i}@f$).
639  xt::xtensor<double, 2> _coeffs;
640 
641  // Number of dofs associated with each cell (sub-)entity
642  //
643  // The dofs of an element are associated with entities of different
644  // topological dimension (vertices, edges, faces, cells). The dofs are
645  // listed in this order, with vertex dofs first. Each entry is the dof
646  // count on the associated entity, as listed by cell::topology.
647  std::vector<std::vector<int>> _num_edofs;
648 
649  // Number of dofs associated with the closure of each cell (sub-)entity
650  std::vector<std::vector<int>> _num_e_closure_dofs;
651 
652  // Dofs associated with each cell (sub-)entity
653  std::vector<std::vector<std::set<int>>> _edofs;
654 
655  // Dofs associated with each cell (sub-)entity
656  std::vector<std::vector<std::set<int>>> _e_closure_dofs;
657 
658  // Entity transformations
659  std::map<cell::type, xt::xtensor<double, 3>> _entity_transformations;
660 
661  // Set of points used for point evaluation
662  // Experimental - currently used for an implementation of
663  // "tabulate_dof_coordinates" Most useful for Lagrange. This may change or go
664  // away. For non-Lagrange elements, these points will be used in combination
665  // with _interpolation_matrix to perform interpolation
666  xt::xtensor<double, 2> _points;
667 
668  // Interpolation points on the cell. The shape is (entity_dim, num
669  // entities of given dimension, num_points, tdim)
670  std::array<std::vector<xt::xtensor<double, 2>>, 4> _x;
671 
673  xt::xtensor<double, 2> _matM;
674 
676  std::array<std::vector<xt::xtensor<double, 3>>, 4> _matM_new;
677 
679  bool _dof_transformations_are_permutations;
680 
682  bool _dof_transformations_are_identity;
683 
687  std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm;
688 
692  std::map<cell::type, std::vector<std::vector<std::size_t>>> _eperm_rev;
693 
695  std::map<cell::type,
696  std::vector<std::tuple<std::vector<std::size_t>, std::vector<double>,
697  xt::xtensor<double, 2>>>>
698  _etrans;
699 
701  std::map<cell::type,
702  std::vector<std::tuple<std::vector<std::size_t>, std::vector<double>,
703  xt::xtensor<double, 2>>>>
704  _etransT;
705 
707  std::map<cell::type,
708  std::vector<std::tuple<std::vector<std::size_t>, std::vector<double>,
709  xt::xtensor<double, 2>>>>
710  _etrans_inv;
711 
713  std::map<cell::type,
714  std::vector<std::tuple<std::vector<std::size_t>, std::vector<double>,
715  xt::xtensor<double, 2>>>>
716  _etrans_invT;
717 };
718 
720 FiniteElement create_element(std::string family, std::string cell, int degree);
721 
723 FiniteElement create_element(element::family family, cell::type cell,
724  int degree);
725 
728 std::string version();
729 
730 //-----------------------------------------------------------------------------
731 template <typename T>
732 void FiniteElement::apply_dof_transformation(const xtl::span<T>& data,
733  int block_size,
734  std::uint32_t cell_info) const
735 {
736  if (_dof_transformations_are_identity)
737  return;
738 
739  if (_cell_tdim >= 2)
740  {
741  // This assumes 3 bits are used per face. This will need updating if
742  // 3D cells with faces with more than 4 sides are implemented
743  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
744  int dofstart
745  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
746 
747  // Transform DOFs on edges
748  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
749  {
750  // Reverse an edge
751  if (cell_info >> (face_start + e) & 1)
752  precompute::apply_matrix(_etrans.at(cell::type::interval)[0], data,
753  dofstart, block_size);
754  dofstart += _num_edofs[1][e];
755  }
756 
757  if (_cell_tdim == 3)
758  {
759  // Permute DOFs on faces
760  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
761  {
762  // Reflect a face
763  if (cell_info >> (3 * f) & 1)
764  precompute::apply_matrix(_etrans.at(_cell_subentity_types[2][f])[1],
765  data, dofstart, block_size);
766 
767  // Rotate a face
768  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
769  precompute::apply_matrix(_etrans.at(_cell_subentity_types[2][f])[0],
770  data, dofstart, block_size);
771  dofstart += _num_edofs[2][f];
772  }
773  }
774  }
775 }
776 //-----------------------------------------------------------------------------
777 template <typename T>
779  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
780 {
781  if (_dof_transformations_are_identity)
782  return;
783 
784  if (_cell_tdim >= 2)
785  {
786  // This assumes 3 bits are used per face. This will need updating if
787  // 3D cells with faces with more than 4 sides are implemented
788  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
789  int dofstart
790  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
791 
792  // Transform DOFs on edges
793  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
794  {
795  // Reverse an edge
796  if (cell_info >> (face_start + e) & 1)
797  precompute::apply_matrix(_etransT.at(cell::type::interval)[0], data,
798  dofstart, block_size);
799  dofstart += _num_edofs[1][e];
800  }
801 
802  if (_cell_tdim == 3)
803  {
804  // Permute DOFs on faces
805  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
806  {
807  // Rotate a face
808  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
809  precompute::apply_matrix(_etransT.at(_cell_subentity_types[2][f])[0],
810  data, dofstart, block_size);
811  // Reflect a face
812  if (cell_info >> (3 * f) & 1)
813  precompute::apply_matrix(_etransT.at(_cell_subentity_types[2][f])[1],
814  data, dofstart, block_size);
815  dofstart += _num_edofs[2][f];
816  }
817  }
818  }
819 }
820 //-----------------------------------------------------------------------------
821 template <typename T>
823  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
824 {
825  if (_dof_transformations_are_identity)
826  return;
827 
828  if (_cell_tdim >= 2)
829  {
830  // This assumes 3 bits are used per face. This will need updating if 3D
831  // cells with faces with more than 4 sides are implemented
832  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
833  int dofstart
834  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
835 
836  // Transform DOFs on edges
837  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
838  {
839  // Reverse an edge
840  if (cell_info >> (face_start + e) & 1)
841  precompute::apply_matrix(_etrans_invT.at(cell::type::interval)[0], data,
842  dofstart, block_size);
843  dofstart += _num_edofs[1][e];
844  }
845 
846  if (_cell_tdim == 3)
847  {
848  // Permute DOFs on faces
849  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
850  {
851  // Reflect a face
852  if (cell_info >> (3 * f) & 1)
854  _etrans_invT.at(_cell_subentity_types[2][f])[1], data, dofstart,
855  block_size);
856 
857  // Rotate a face
858  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
860  _etrans_invT.at(_cell_subentity_types[2][f])[0], data, dofstart,
861  block_size);
862  dofstart += _num_edofs[2][f];
863  }
864  }
865  }
866 }
867 //-----------------------------------------------------------------------------
868 template <typename T>
870  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
871 {
872  if (_dof_transformations_are_identity)
873  return;
874 
875  if (_cell_tdim >= 2)
876  {
877  // This assumes 3 bits are used per face. This will need updating if 3D
878  // cells with faces with more than 4 sides are implemented
879  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
880  int dofstart
881  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
882 
883  // Transform DOFs on edges
884  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
885  {
886  // Reverse an edge
887  if (cell_info >> (face_start + e) & 1)
888  precompute::apply_matrix(_etrans_inv.at(cell::type::interval)[0], data,
889  dofstart, block_size);
890  dofstart += _num_edofs[1][e];
891  }
892 
893  if (_cell_tdim == 3)
894  {
895  // Permute DOFs on faces
896  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
897  {
898  // Rotate a face
899  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
901  _etrans_inv.at(_cell_subentity_types[2][f])[0], data, dofstart,
902  block_size);
903  // Reflect a face
904  if (cell_info >> (3 * f) & 1)
906  _etrans_inv.at(_cell_subentity_types[2][f])[1], data, dofstart,
907  block_size);
908  dofstart += _num_edofs[2][f];
909  }
910  }
911  }
912 }
913 //-----------------------------------------------------------------------------
914 template <typename T>
916  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
917 {
918  if (_dof_transformations_are_identity)
919  return;
920 
921  if (_cell_tdim >= 2)
922  {
923  // This assumes 3 bits are used per face. This will need updating if
924  // 3D cells with faces with more than 4 sides are implemented
925  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
926  int dofstart
927  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
928 
929  // Transform DOFs on edges
930  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
931  {
932  // Reverse an edge
933  if (cell_info >> (face_start + e) & 1)
935  _etrans.at(cell::type::interval)[0], data, dofstart, block_size);
936  dofstart += _num_edofs[1][e];
937  }
938 
939  if (_cell_tdim == 3)
940  {
941  // Permute DOFs on faces
942  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
943  {
944  // Reflect a face
945  if (cell_info >> (3 * f) & 1)
947  _etrans.at(_cell_subentity_types[2][f])[1], data, dofstart,
948  block_size);
949 
950  // Rotate a face
951  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
953  _etrans.at(_cell_subentity_types[2][f])[0], data, dofstart,
954  block_size);
955  dofstart += _num_edofs[2][f];
956  }
957  }
958  }
959 }
960 //-----------------------------------------------------------------------------
961 template <typename T>
963  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
964 {
965  if (_dof_transformations_are_identity)
966  return;
967 
968  if (_cell_tdim >= 2)
969  {
970  // This assumes 3 bits are used per face. This will need updating if
971  // 3D cells with faces with more than 4 sides are implemented
972  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
973  int dofstart
974  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
975 
976  // Transform DOFs on edges
977  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
978  {
979  // Reverse an edge
980  if (cell_info >> (face_start + e) & 1)
982  _etrans_invT.at(cell::type::interval)[0], data, dofstart,
983  block_size);
984  dofstart += _num_edofs[1][e];
985  }
986 
987  if (_cell_tdim == 3)
988  {
989  // Permute DOFs on faces
990  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
991  {
992  // Reflect a face
993  if (cell_info >> (3 * f) & 1)
995  _etrans_invT.at(_cell_subentity_types[2][f])[1], data, dofstart,
996  block_size);
997 
998  // Rotate a face
999  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1001  _etrans_invT.at(_cell_subentity_types[2][f])[0], data, dofstart,
1002  block_size);
1003  dofstart += _num_edofs[2][f];
1004  }
1005  }
1006  }
1007 }
1008 //-----------------------------------------------------------------------------
1009 template <typename T>
1011  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1012 {
1013  if (_dof_transformations_are_identity)
1014  return;
1015 
1016  if (_cell_tdim >= 2)
1017  {
1018  // This assumes 3 bits are used per face. This will need updating if
1019  // 3D cells with faces with more than 4 sides are implemented
1020  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1021  int dofstart
1022  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1023 
1024  // Transform DOFs on edges
1025  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1026  {
1027  // Reverse an edge
1028  if (cell_info >> (face_start + e) & 1)
1030  _etransT.at(cell::type::interval)[0], data, dofstart, block_size);
1031  dofstart += _num_edofs[1][e];
1032  }
1033 
1034  if (_cell_tdim == 3)
1035  {
1036  // Permute DOFs on faces
1037  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1038  {
1039  // Rotate a face
1040  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1042  _etransT.at(_cell_subentity_types[2][f])[0], data, dofstart,
1043  block_size);
1044 
1045  // Reflect a face
1046  if (cell_info >> (3 * f) & 1)
1048  _etransT.at(_cell_subentity_types[2][f])[1], data, dofstart,
1049  block_size);
1050 
1051  dofstart += _num_edofs[2][f];
1052  }
1053  }
1054  }
1055 }
1056 //-----------------------------------------------------------------------------
1057 template <typename T>
1059  const xtl::span<T>& data, int block_size, std::uint32_t cell_info) const
1060 {
1061  if (_dof_transformations_are_identity)
1062  return;
1063 
1064  if (_cell_tdim >= 2)
1065  {
1066  // This assumes 3 bits are used per face. This will need updating if
1067  // 3D cells with faces with more than 4 sides are implemented
1068  int face_start = _cell_tdim == 3 ? 3 * _num_edofs[2].size() : 0;
1069  int dofstart
1070  = std::accumulate(_num_edofs[0].cbegin(), _num_edofs[0].cend(), 0);
1071 
1072  // Transform DOFs on edges
1073  for (std::size_t e = 0; e < _num_edofs[1].size(); ++e)
1074  {
1075  // Reverse an edge
1076  if (cell_info >> (face_start + e) & 1)
1078  _etrans_inv.at(cell::type::interval)[0], data, dofstart,
1079  block_size);
1080  dofstart += _num_edofs[1][e];
1081  }
1082 
1083  if (_cell_tdim == 3)
1084  {
1085  // Permute DOFs on faces
1086  for (std::size_t f = 0; f < _num_edofs[2].size(); ++f)
1087  {
1088  // Rotate a face
1089  for (std::uint32_t r = 0; r < (cell_info >> (3 * f + 1) & 3); ++r)
1091  _etrans_inv.at(_cell_subentity_types[2][f])[0], data, dofstart,
1092  block_size);
1093 
1094  // Reflect a face
1095  if (cell_info >> (3 * f) & 1)
1097  _etrans_inv.at(_cell_subentity_types[2][f])[1], data, dofstart,
1098  block_size);
1099 
1100  dofstart += _num_edofs[2][f];
1101  }
1102  }
1103  }
1104 }
1105 //-----------------------------------------------------------------------------
1106 
1107 } // namespace basix
basix::FiniteElement::permute_dofs
void permute_dofs(const xtl::span< std::int32_t > &dofs, std::uint32_t cell_info) const
Definition: finite-element.cpp:671
basix::FiniteElement::entity_dofs
const std::vector< std::vector< std::set< int > > > & entity_dofs() const
Definition: finite-element.cpp:523
basix::FiniteElement::dof_transformations_are_identity
bool dof_transformations_are_identity() const
Definition: finite-element.cpp:507
basix::FiniteElement::tabulate
xt::xtensor< double, 4 > tabulate(int nd, const xt::xarray< double > &x) const
Definition: finite-element.cpp:541
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:778
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:732
basix::FiniteElement::degree
int degree() const
Definition: finite-element.cpp:483
basix::FiniteElement::~FiniteElement
~FiniteElement()=default
Destructor.
basix::FiniteElement::num_entity_dofs
const std::vector< std::vector< int > > & num_entity_dofs() const
Definition: finite-element.cpp:517
basix::cell::type
type
Cell type.
Definition: cell.h:16
basix::FiniteElement::points
const xt::xtensor< double, 2 > & points() const
Definition: finite-element.cpp:648
basix::FiniteElement::entity_transformations
std::map< cell::type, xt::xtensor< double, 3 > > entity_transformations() const
Return the entity dof transformation matricess.
Definition: finite-element.cpp:773
basix::FiniteElement::cell_type
cell::type cell_type() const
Definition: finite-element.cpp:481
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:962
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:869
basix::FiniteElement::mapping_type
maps::type mapping_type() const
Definition: finite-element.cpp:500
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:1010
basix::FiniteElement::dof_transformations_are_permutations
bool dof_transformations_are_permutations() const
Definition: finite-element.cpp:502
basix::compute_expansion_coefficients
xt::xtensor< double, 3 > compute_expansion_coefficients(cell::type cell_type, const xt::xtensor< double, 2 > &B, const std::vector< std::vector< xt::xtensor< double, 3 >>> &M, const std::vector< std::vector< xt::xtensor< double, 2 >>> &x, int degree, double kappa_tol=0.0)
Definition: finite-element.cpp:152
basix::FiniteElement::family
element::family family() const
Definition: finite-element.cpp:498
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:1058
basix::FiniteElement::dim
int dim() const
Definition: finite-element.cpp:496
basix::FiniteElement::map_pull_back
xt::xtensor< double, 3 > map_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:661
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:915
basix::maps::apply_map
void apply_map(O &&u, const P &U, const Mat0 &J, double detJ, const Mat1 &K, maps::type map_type)
Definition: maps.h:130
basix::maps::type
type
Cell type.
Definition: maps.h:19
basix::FiniteElement::value_shape
const std::vector< int > & value_shape() const
Definition: finite-element.cpp:491
basix::FiniteElement::entity_closure_dofs
const std::vector< std::vector< std::set< int > > > & entity_closure_dofs() const
Definition: finite-element.cpp:535
basix::FiniteElement::base_transformations
xt::xtensor< double, 3 > base_transformations() const
Definition: finite-element.cpp:588
basix::FiniteElement::num_entity_closure_dofs
const std::vector< std::vector< int > > & num_entity_closure_dofs() const
Definition: finite-element.cpp:529
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:340
basix::FiniteElement
Definition: finite-element.h:193
basix::FiniteElement::interpolation_matrix
const xt::xtensor< double, 2 > & interpolation_matrix() const
Definition: finite-element.cpp:512
basix::FiniteElement::map_push_forward
xt::xtensor< double, 3 > map_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:650
basix::version
std::string version()
Definition: finite-element.cpp:778
basix::FiniteElement::value_size
int value_size() const
Definition: finite-element.cpp:485
basix::FiniteElement::map_type
maps::type map_type
Element map type.
Definition: finite-element.h:610
basix::FiniteElement::map_push_forward_m
void map_push_forward_m(const O &U, const P &J, const Q &detJ, const S &K, T &&u) const
Definition: finite-element.h:335
basix
Placeholder.
Definition: brezzi-douglas-marini.h:9
basix::FiniteElement::unpermute_dofs
void unpermute_dofs(const xtl::span< std::int32_t > &dofs, std::uint32_t cell_info) const
Definition: finite-element.cpp:722
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:311
basix::FiniteElement::num_points
int num_points() const
Return the number of interpolation points.
Definition: finite-element.cpp:646
basix::FiniteElement::operator=
FiniteElement & operator=(const FiniteElement &element)=default
Assignment operator.
basix::FiniteElement::map_pull_back_m
void map_pull_back_m(const O &u, const P &J, const Q &detJ, const S &K, T &&U) const
Definition: finite-element.h:379
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:822
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, 3 > &coeffs, const std::map< cell::type, xt::xtensor< double, 3 >> &entity_transformations, 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=maps::type::identity)
Definition: finite-element.cpp:252
basix::create_element
FiniteElement create_element(std::string family, std::string cell, int degree)
Create an element by name.
Definition: finite-element.cpp:77