Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
reference_cell.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2020 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_tria_reference_cell_h
17#define dealii_tria_reference_cell_h
18
19#include <deal.II/base/config.h>
20
24#include <deal.II/base/tensor.h>
26
27
29
30// Forward declarations
31#ifndef DOXYGEN
32template <int dim, int spacedim>
33class Mapping;
34
35template <int dim>
36class Quadrature;
37
38class ReferenceCell;
39#endif
40
41
42namespace internal
43{
44 namespace ReferenceCell
45 {
59 DEAL_II_CONSTEXPR ::ReferenceCell
60 make_reference_cell_from_int(const std::uint8_t kind);
61
62 } // namespace ReferenceCell
63} // namespace internal
64
65
66
88{
89public:
97 static ReferenceCell
98 n_vertices_to_type(const int dim, const unsigned int n_vertices);
99
110
121 bool
122 is_hyper_cube() const;
123
127 bool
128 is_simplex() const;
129
134 unsigned int
135 get_dimension() const;
136
150 template <int dim>
151 double
152 d_linear_shape_function(const Point<dim> &xi, const unsigned int i) const;
153
158 template <int dim>
161 const unsigned int i) const;
162
172 template <int dim, int spacedim = dim>
173 std::unique_ptr<Mapping<dim, spacedim>>
174 get_default_mapping(const unsigned int degree) const;
175
187 template <int dim, int spacedim = dim>
190
199 template <int dim>
201 get_gauss_type_quadrature(const unsigned n_points_1D) const;
202
209 template <int dim>
210 const Quadrature<dim> &
212
227 unsigned int
228 n_vertices() const;
229
235 vertex_indices() const;
236
242 unsigned int
243 n_lines() const;
244
250 line_indices() const;
251
257 unsigned int
258 n_faces() const;
259
265 face_indices() const;
266
280 face_reference_cell(const unsigned int face_no) const;
281
332 unsigned int
333 child_cell_on_face(const unsigned int face_n,
334 const unsigned int subface_n,
335 const unsigned char face_orientation = 1) const;
336
346 std::array<unsigned int, 2>
347 standard_vertex_to_face_and_vertex_index(const unsigned int vertex) const;
348
358 std::array<unsigned int, 2>
359 standard_line_to_face_and_line_index(const unsigned int line) const;
360
364 unsigned int
365 face_to_cell_lines(const unsigned int face,
366 const unsigned int line,
367 const unsigned char face_orientation) const;
368
372 unsigned int
373 face_to_cell_vertices(const unsigned int face,
374 const unsigned int vertex,
375 const unsigned char face_orientation) const;
376
380 unsigned int
381 standard_to_real_face_vertex(const unsigned int vertex,
382 const unsigned int face,
383 const unsigned char face_orientation) const;
384
388 unsigned int
389 standard_to_real_face_line(const unsigned int line,
390 const unsigned int face,
391 const unsigned char face_orientation) const;
392
403 bool
404 standard_vs_true_line_orientation(const unsigned int line,
405 const unsigned char face_orientation,
406 const unsigned char line_orientation) const;
407
418 /*
419 * Return @f$i@f$-th unit tangential vector of a face of the reference cell.
420 * The vectors are arranged such that the
421 * cross product between the two vectors returns the unit normal vector.
422 *
423 * @pre @f$i@f$ must be between zero and `dim-1`.
424 */
425 template <int dim>
427 unit_tangential_vectors(const unsigned int face_no,
428 const unsigned int i) const;
429
433 template <int dim>
435 unit_normal_vectors(const unsigned int face_no) const;
436
441 template <typename T, std::size_t N>
442 unsigned char
443 compute_orientation(const std::array<T, N> &vertices_0,
444 const std::array<T, N> &vertices_1) const;
445
449 template <typename T, std::size_t N>
450 std::array<T, N>
451 permute_according_orientation(const std::array<T, N> &vertices,
452 const unsigned int orientation) const;
453
458 faces_for_given_vertex(const unsigned int vertex_index) const;
459
472 unsigned int
473 exodusii_vertex_to_deal_vertex(const unsigned int vertex_n) const;
474
478 unsigned int
479 exodusii_face_to_deal_face(const unsigned int face_n) const;
480
484 unsigned int
485 unv_vertex_to_deal_vertex(const unsigned int vertex_n) const;
486
490 unsigned int
491 vtk_linear_type() const;
492
497 unsigned int
498 vtk_quadratic_type() const;
499
504 unsigned int
505 vtk_lagrange_type() const;
506
520 std::string
521 to_string() const;
522
526 constexpr operator std::uint8_t() const;
527
531 constexpr bool
532 operator==(const ReferenceCell &type) const;
533
537 constexpr bool
538 operator!=(const ReferenceCell &type) const;
539
545 template <class Archive>
546 void
547 serialize(Archive &archive, const unsigned int /*version*/);
548
553private:
557 std::uint8_t kind;
558
565 constexpr ReferenceCell(const std::uint8_t kind);
566
573};
574
575
576
577inline constexpr ReferenceCell::ReferenceCell(const std::uint8_t kind)
578 : kind(kind)
579{}
580
581
582
583inline constexpr ReferenceCell::operator std::uint8_t() const
584{
585 return kind;
586}
587
588
589
590inline constexpr bool
592{
593 return kind == type.kind;
594}
595
596
597
598inline constexpr bool
600{
601 return kind != type.kind;
602}
603
604
605
606namespace internal
607{
608 namespace ReferenceCell
609 {
610 inline DEAL_II_CONSTEXPR ::ReferenceCell
611 make_reference_cell_from_int(const std::uint8_t kind)
612 {
613 // Make sure these are the only indices from which objects can be
614 // created.
615 Assert((kind == static_cast<std::uint8_t>(-1)) || (kind < 8),
617
618 // Call the private constructor, which we can from here because this
619 // function is a 'friend'.
620 return {kind};
621 }
622 } // namespace ReferenceCell
623} // namespace internal
624
625
626
635{
654 static_cast<std::uint8_t>(-1));
655
661 template <int dim>
662 constexpr const ReferenceCell &
663 get_simplex();
664
670 template <int dim>
671 constexpr const ReferenceCell &
673} // namespace ReferenceCells
674
675
676
680{}
681
682
683
684template <class Archive>
685inline void
686ReferenceCell::serialize(Archive &archive, const unsigned int /*version*/)
687{
688 archive &kind;
689}
690
691
692
694ReferenceCell::faces_for_given_vertex(const unsigned int vertex) const
695{
696 if (*this == ReferenceCells::Line)
697 {
699 return {&GeometryInfo<2>::vertex_to_face[vertex][0], 1};
700 }
701 else if (*this == ReferenceCells::Quadrilateral)
702 {
704 return {&GeometryInfo<2>::vertex_to_face[vertex][0], 2};
705 }
706 else if (*this == ReferenceCells::Hexahedron)
707 {
709 return {&GeometryInfo<3>::vertex_to_face[vertex][0], 3};
710 }
711 else if (*this == ReferenceCells::Triangle)
712 {
713 AssertIndexRange(vertex, 3);
714 static const ndarray<unsigned int, 3, 2> table = {
715 {{{0, 2}}, {{0, 1}}, {{1, 2}}}};
716
717 return table[vertex];
718 }
719 else if (*this == ReferenceCells::Tetrahedron)
720 {
721 AssertIndexRange(vertex, 4);
722 static const ndarray<unsigned int, 4, 3> table = {
723 {{{0, 1, 2}}, {{0, 1, 3}}, {{0, 2, 3}}, {{1, 2, 3}}}};
724
725 return table[vertex];
726 }
727 else if (*this == ReferenceCells::Wedge)
728 {
729 AssertIndexRange(vertex, 6);
730 static const ndarray<unsigned int, 6, 3> table = {{{{0, 2, 4}},
731 {{0, 2, 3}},
732 {{0, 3, 4}},
733 {{1, 2, 4}},
734 {{1, 2, 3}},
735 {{1, 3, 4}}}};
736
737 return table[vertex];
738 }
739 else if (*this == ReferenceCells::Pyramid)
740 {
741 AssertIndexRange(vertex, 5);
742 static const unsigned int X = numbers::invalid_unsigned_int;
743 static const ndarray<unsigned int, 5, 4> table = {{{{0, 1, 3, X}},
744 {{0, 2, 3, X}},
745 {{0, 1, 4, X}},
746 {{0, 2, 4, X}},
747 {{1, 2, 3, 4}}}};
748
749 return {&table[vertex][0], vertex == 4 ? 4u : 3u};
750 }
751
752 Assert(false, ExcNotImplemented());
753
754 return {};
755}
756
757
758
759inline bool
761{
762 return (*this == ReferenceCells::Vertex || *this == ReferenceCells::Line ||
765}
766
767
768
769inline bool
771{
772 return (*this == ReferenceCells::Vertex || *this == ReferenceCells::Line ||
773 *this == ReferenceCells::Triangle ||
775}
776
777
778
779inline unsigned int
781{
782 if (*this == ReferenceCells::Vertex)
783 return 0;
784 else if (*this == ReferenceCells::Line)
785 return 1;
786 else if ((*this == ReferenceCells::Triangle) ||
788 return 2;
789 else if ((*this == ReferenceCells::Tetrahedron) ||
790 (*this == ReferenceCells::Pyramid) ||
791 (*this == ReferenceCells::Wedge) ||
793 return 3;
794
795 Assert(false, ExcNotImplemented());
797}
798
799
800
801inline unsigned int
803{
804 if (*this == ReferenceCells::Vertex)
805 return 1;
806 else if (*this == ReferenceCells::Line)
807 return 2;
808 else if (*this == ReferenceCells::Triangle)
809 return 3;
810 else if (*this == ReferenceCells::Quadrilateral)
811 return 4;
812 else if (*this == ReferenceCells::Tetrahedron)
813 return 4;
814 else if (*this == ReferenceCells::Pyramid)
815 return 5;
816 else if (*this == ReferenceCells::Wedge)
817 return 6;
818 else if (*this == ReferenceCells::Hexahedron)
819 return 8;
820
821 Assert(false, ExcNotImplemented());
823}
824
825
826
827inline unsigned int
829{
830 if (*this == ReferenceCells::Vertex)
831 return 0;
832 else if (*this == ReferenceCells::Line)
833 return 1;
834 else if (*this == ReferenceCells::Triangle)
835 return 3;
836 else if (*this == ReferenceCells::Quadrilateral)
837 return 4;
838 else if (*this == ReferenceCells::Tetrahedron)
839 return 6;
840 else if (*this == ReferenceCells::Pyramid)
841 return 7;
842 else if (*this == ReferenceCells::Wedge)
843 return 9;
844 else if (*this == ReferenceCells::Hexahedron)
845 return 12;
846
847 Assert(false, ExcNotImplemented());
849}
850
851
852
853inline unsigned int
855{
856 if (*this == ReferenceCells::Vertex)
857 return 0;
858 else if (*this == ReferenceCells::Line)
859 return 2;
860 else if (*this == ReferenceCells::Triangle)
861 return 3;
862 else if (*this == ReferenceCells::Quadrilateral)
863 return 4;
864 else if (*this == ReferenceCells::Tetrahedron)
865 return 4;
866 else if (*this == ReferenceCells::Pyramid)
867 return 5;
868 else if (*this == ReferenceCells::Wedge)
869 return 5;
870 else if (*this == ReferenceCells::Hexahedron)
871 return 6;
872
873 Assert(false, ExcNotImplemented());
875}
876
877
878
881{
882 return {0U, n_vertices()};
883}
884
885
886
889{
890 return {0U, n_lines()};
891}
892
893
894
897{
898 return {0U, n_faces()};
899}
900
901
902
903inline ReferenceCell
904ReferenceCell::face_reference_cell(const unsigned int face_no) const
905{
906 AssertIndexRange(face_no, n_faces());
907
908 if (*this == ReferenceCells::Vertex)
910 else if (*this == ReferenceCells::Line)
912 else if (*this == ReferenceCells::Triangle)
914 else if (*this == ReferenceCells::Quadrilateral)
916 else if (*this == ReferenceCells::Tetrahedron)
918 else if (*this == ReferenceCells::Pyramid)
919 {
920 if (face_no == 0)
922 else
924 }
925 else if (*this == ReferenceCells::Wedge)
926 {
927 if (face_no > 1)
929 else
931 }
932 else if (*this == ReferenceCells::Hexahedron)
934
935 Assert(false, ExcNotImplemented());
937}
938
939
940
941inline unsigned int
943 const unsigned int face,
944 const unsigned int subface,
945 const unsigned char face_orientation_raw) const
946{
947 AssertIndexRange(face, n_faces());
948
949 if (*this == ReferenceCells::Vertex)
950 {
951 Assert(false, ExcNotImplemented());
952 }
953 else if (*this == ReferenceCells::Line)
954 {
955 Assert(false, ExcNotImplemented());
956 }
957 else if (*this == ReferenceCells::Triangle)
958 {
959 static const ndarray<unsigned int, 3, 2> subcells = {
960 {{{0, 1}}, {{1, 2}}, {{2, 0}}}};
961
962 return subcells[face][subface];
963 }
964 else if (*this == ReferenceCells::Quadrilateral)
965 {
966 const bool face_orientation = Utilities::get_bit(face_orientation_raw, 0);
967 const bool face_flip = Utilities::get_bit(face_orientation_raw, 2);
968 const bool face_rotation = Utilities::get_bit(face_orientation_raw, 1);
969
972 face,
973 subface,
974 face_orientation,
975 face_flip,
976 face_rotation);
977 }
978 else if (*this == ReferenceCells::Tetrahedron)
979 {
980 Assert(false, ExcNotImplemented());
981 }
982 else if (*this == ReferenceCells::Pyramid)
983 {
984 Assert(false, ExcNotImplemented());
985 }
986 else if (*this == ReferenceCells::Wedge)
987 {
988 Assert(false, ExcNotImplemented());
989 }
990 else if (*this == ReferenceCells::Hexahedron)
991 {
992 const bool face_orientation = Utilities::get_bit(face_orientation_raw, 0);
993 const bool face_flip = Utilities::get_bit(face_orientation_raw, 2);
994 const bool face_rotation = Utilities::get_bit(face_orientation_raw, 1);
995
998 face,
999 subface,
1000 face_orientation,
1001 face_flip,
1002 face_rotation);
1003 }
1004
1005 Assert(false, ExcNotImplemented());
1006 return {};
1007}
1008
1009
1010
1011inline std::array<unsigned int, 2>
1013 const unsigned int vertex) const
1014{
1015 AssertIndexRange(vertex, n_vertices());
1016
1017 if (*this == ReferenceCells::Vertex)
1018 {
1019 Assert(false, ExcNotImplemented());
1020 }
1021 else if (*this == ReferenceCells::Line)
1022 {
1023 Assert(false, ExcNotImplemented());
1024 }
1025 else if (*this == ReferenceCells::Triangle)
1026 {
1027 static const ndarray<unsigned int, 3, 2> table = {
1028 {{{0, 0}}, {{0, 1}}, {{1, 1}}}};
1029
1030 return table[vertex];
1031 }
1032 else if (*this == ReferenceCells::Quadrilateral)
1033 {
1035 }
1036 else if (*this == ReferenceCells::Tetrahedron)
1037 {
1038 static const ndarray<unsigned int, 4, 2> table = {
1039 {{{0, 0}}, {{0, 1}}, {{0, 2}}, {{1, 2}}}};
1040
1041 return table[vertex];
1042 }
1043 else if (*this == ReferenceCells::Pyramid)
1044 {
1045 static const ndarray<unsigned int, 5, 2> table = {
1046 {{{0, 0}}, {{0, 1}}, {{0, 2}}, {{0, 3}}, {{1, 2}}}};
1047
1048 return table[vertex];
1049 }
1050 else if (*this == ReferenceCells::Wedge)
1051 {
1052 static const ndarray<unsigned int, 6, 2> table = {
1053 {{{0, 1}}, {{0, 0}}, {{0, 2}}, {{1, 0}}, {{1, 1}}, {{1, 2}}}};
1054
1055 return table[vertex];
1056 }
1057 else if (*this == ReferenceCells::Hexahedron)
1058 {
1060 }
1061
1062 Assert(false, ExcNotImplemented());
1063 return {};
1064}
1065
1066
1067
1068inline std::array<unsigned int, 2>
1070 const unsigned int line) const
1071{
1072 AssertIndexRange(line, n_lines());
1073
1074 if (*this == ReferenceCells::Vertex)
1075 {
1076 Assert(false, ExcNotImplemented());
1077 }
1078 else if (*this == ReferenceCells::Line)
1079 {
1080 Assert(false, ExcNotImplemented());
1081 }
1082 else if (*this == ReferenceCells::Triangle)
1083 {
1084 Assert(false, ExcNotImplemented());
1085 }
1086 else if (*this == ReferenceCells::Quadrilateral)
1087 {
1088 Assert(false, ExcNotImplemented());
1089 }
1090 else if (*this == ReferenceCells::Tetrahedron)
1091 {
1092 static const std::array<unsigned int, 2> table[6] = {
1093 {{0, 0}}, {{0, 1}}, {{0, 2}}, {{1, 1}}, {{1, 2}}, {{2, 1}}};
1094
1095 return table[line];
1096 }
1097 else if (*this == ReferenceCells::Pyramid)
1098 {
1099 static const std::array<unsigned int, 2> table[8] = {{{0, 0}},
1100 {{0, 1}},
1101 {{0, 2}},
1102 {{0, 3}},
1103 {{1, 2}},
1104 {{2, 1}},
1105 {{1, 1}},
1106 {{2, 2}}};
1107
1108 return table[line];
1109 }
1110 else if (*this == ReferenceCells::Wedge)
1111 {
1112 static const std::array<unsigned int, 2> table[9] = {{{0, 0}},
1113 {{0, 2}},
1114 {{0, 1}},
1115 {{1, 0}},
1116 {{1, 1}},
1117 {{1, 2}},
1118 {{2, 0}},
1119 {{2, 1}},
1120 {{3, 1}}};
1121
1122 return table[line];
1123 }
1124 else if (*this == ReferenceCells::Hexahedron)
1125 {
1127 }
1128
1129 Assert(false, ExcNotImplemented());
1130 return {};
1131}
1132
1133
1134
1135inline unsigned int
1137 const unsigned int line,
1138 const unsigned char face_orientation) const
1139{
1140 AssertIndexRange(face, n_faces());
1142
1143 if (*this == ReferenceCells::Vertex)
1144 {
1145 Assert(false, ExcNotImplemented());
1146 }
1147 else if (*this == ReferenceCells::Line)
1148 {
1150 face,
1151 line,
1152 Utilities::get_bit(face_orientation, 0),
1153 Utilities::get_bit(face_orientation, 2),
1154 Utilities::get_bit(face_orientation, 1));
1155 }
1156 else if (*this == ReferenceCells::Triangle)
1157 {
1158 return face;
1159 }
1160 else if (*this == ReferenceCells::Quadrilateral)
1161 {
1163 face,
1164 line,
1165 Utilities::get_bit(face_orientation, 0),
1166 Utilities::get_bit(face_orientation, 2),
1167 Utilities::get_bit(face_orientation, 1));
1168 }
1169 else if (*this == ReferenceCells::Tetrahedron)
1170 {
1171 const static ndarray<unsigned int, 4, 3> table = {
1172 {{{0, 1, 2}}, {{0, 3, 4}}, {{2, 5, 3}}, {{1, 4, 5}}}};
1173
1174 return table[face]
1175 [standard_to_real_face_line(line, face, face_orientation)];
1176 }
1177 else if (*this == ReferenceCells::Pyramid)
1178 {
1179 Assert(false, ExcNotImplemented());
1180 }
1181 else if (*this == ReferenceCells::Wedge)
1182 {
1183 Assert(false, ExcNotImplemented());
1184 }
1185 else if (*this == ReferenceCells::Hexahedron)
1186 {
1188 face,
1189 line,
1190 Utilities::get_bit(face_orientation, 0),
1191 Utilities::get_bit(face_orientation, 2),
1192 Utilities::get_bit(face_orientation, 1));
1193 }
1194
1195 Assert(false, ExcNotImplemented());
1197}
1198
1199
1200
1201inline unsigned int
1203 const unsigned int vertex,
1204 const unsigned char face_orientation) const
1205{
1206 AssertIndexRange(face, n_faces());
1208
1209 if (*this == ReferenceCells::Vertex)
1210 {
1211 Assert(false, ExcNotImplemented());
1212 }
1213 else if (*this == ReferenceCells::Line)
1214 {
1216 face,
1217 vertex,
1218 Utilities::get_bit(face_orientation, 0),
1219 Utilities::get_bit(face_orientation, 2),
1220 Utilities::get_bit(face_orientation, 1));
1221 }
1222 else if (*this == ReferenceCells::Triangle)
1223 {
1224 static const ndarray<unsigned int, 3, 2> table = {
1225 {{{0, 1}}, {{1, 2}}, {{2, 0}}}};
1226
1227 return table[face][face_orientation ? vertex : (1 - vertex)];
1228 }
1229 else if (*this == ReferenceCells::Quadrilateral)
1230 {
1232 face,
1233 vertex,
1234 Utilities::get_bit(face_orientation, 0),
1235 Utilities::get_bit(face_orientation, 2),
1236 Utilities::get_bit(face_orientation, 1));
1237 }
1238 else if (*this == ReferenceCells::Tetrahedron)
1239 {
1240 static const ndarray<unsigned int, 4, 3> table = {
1241 {{{0, 1, 2}}, {{1, 0, 3}}, {{0, 2, 3}}, {{2, 1, 3}}}};
1242
1243 return table[face][standard_to_real_face_vertex(
1244 vertex, face, face_orientation)];
1245 }
1246 else if (*this == ReferenceCells::Pyramid)
1247 {
1248 constexpr auto X = numbers::invalid_unsigned_int;
1249 static const ndarray<unsigned int, 5, 4> table = {{{{0, 1, 2, 3}},
1250 {{0, 2, 4, X}},
1251 {{3, 1, 4, X}},
1252 {{1, 0, 4, X}},
1253 {{2, 3, 4, X}}}};
1254
1255 return table[face][standard_to_real_face_vertex(
1256 vertex, face, face_orientation)];
1257 }
1258 else if (*this == ReferenceCells::Wedge)
1259 {
1260 constexpr auto X = numbers::invalid_unsigned_int;
1261 static const ndarray<unsigned int, 6, 4> table = {{{{1, 0, 2, X}},
1262 {{3, 4, 5, X}},
1263 {{0, 1, 3, 4}},
1264 {{1, 2, 4, 5}},
1265 {{2, 0, 5, 3}}}};
1266
1267 return table[face][standard_to_real_face_vertex(
1268 vertex, face, face_orientation)];
1269 }
1270 else if (*this == ReferenceCells::Hexahedron)
1271 {
1273 face,
1274 vertex,
1275 Utilities::get_bit(face_orientation, 0),
1276 Utilities::get_bit(face_orientation, 2),
1277 Utilities::get_bit(face_orientation, 1));
1278 }
1279
1280 Assert(false, ExcNotImplemented());
1282}
1283
1284
1285
1286inline unsigned int
1288 const unsigned int vertex,
1289 const unsigned int face,
1290 const unsigned char face_orientation) const
1291{
1292 AssertIndexRange(face, n_faces());
1294
1295 if (*this == ReferenceCells::Vertex)
1296 {
1297 Assert(false, ExcNotImplemented());
1298 }
1299 else if (*this == ReferenceCells::Line)
1300 {
1301 Assert(false, ExcNotImplemented());
1302 }
1303 else if (*this == ReferenceCells::Triangle)
1304 {
1305 static const ndarray<unsigned int, 2, 2> table = {{{{1, 0}}, {{0, 1}}}};
1306
1307 return table[face_orientation][vertex];
1308 }
1309 else if (*this == ReferenceCells::Quadrilateral)
1310 {
1312 face_orientation);
1313 }
1314 else if (*this == ReferenceCells::Tetrahedron)
1315 {
1316 static const ndarray<unsigned int, 6, 3> table = {{{{0, 2, 1}},
1317 {{0, 1, 2}},
1318 {{2, 1, 0}},
1319 {{1, 2, 0}},
1320 {{1, 0, 2}},
1321 {{2, 0, 1}}}};
1322
1323 return table[face_orientation][vertex];
1324 }
1325 else if (*this == ReferenceCells::Pyramid)
1326 {
1327 if (face == 0) // The quadrilateral face
1328 {
1330 vertex,
1331 Utilities::get_bit(face_orientation, 0),
1332 Utilities::get_bit(face_orientation, 2),
1333 Utilities::get_bit(face_orientation, 1));
1334 }
1335 else // One of the triangular faces
1336 {
1337 static const ndarray<unsigned int, 6, 3> table = {{{{0, 2, 1}},
1338 {{0, 1, 2}},
1339 {{2, 1, 0}},
1340 {{1, 2, 0}},
1341 {{1, 0, 2}},
1342 {{2, 0, 1}}}};
1343
1344 return table[face_orientation][vertex];
1345 }
1346 }
1347 else if (*this == ReferenceCells::Wedge)
1348 {
1349 if (face > 1) // One of the quadrilateral faces
1350 {
1352 vertex,
1353 Utilities::get_bit(face_orientation, 0),
1354 Utilities::get_bit(face_orientation, 2),
1355 Utilities::get_bit(face_orientation, 1));
1356 }
1357 else // One of the triangular faces
1358 {
1359 static const ndarray<unsigned int, 6, 3> table = {{{{0, 2, 1}},
1360 {{0, 1, 2}},
1361 {{2, 1, 0}},
1362 {{1, 2, 0}},
1363 {{1, 0, 2}},
1364 {{2, 0, 1}}}};
1365
1366 return table[face_orientation][vertex];
1367 }
1368 }
1369 else if (*this == ReferenceCells::Hexahedron)
1370 {
1372 vertex,
1373 Utilities::get_bit(face_orientation, 0),
1374 Utilities::get_bit(face_orientation, 2),
1375 Utilities::get_bit(face_orientation, 1));
1376 }
1377
1378 Assert(false, ExcNotImplemented());
1380}
1381
1382
1383
1384inline unsigned int
1386 const unsigned int line,
1387 const unsigned int face,
1388 const unsigned char face_orientation) const
1389{
1390 AssertIndexRange(face, n_faces());
1392
1393 if (*this == ReferenceCells::Vertex)
1394 {
1395 Assert(false, ExcNotImplemented());
1396 }
1397 else if (*this == ReferenceCells::Line)
1398 {
1399 Assert(false, ExcNotImplemented());
1400 }
1401 else if (*this == ReferenceCells::Triangle)
1402 {
1403 Assert(false, ExcNotImplemented());
1404 }
1405 else if (*this == ReferenceCells::Quadrilateral)
1406 {
1407 Assert(false, ExcNotImplemented());
1408 }
1409 else if (*this == ReferenceCells::Tetrahedron)
1410 {
1411 static const ndarray<unsigned int, 6, 3> table = {{{{2, 1, 0}},
1412 {{0, 1, 2}},
1413 {{1, 0, 2}},
1414 {{1, 2, 0}},
1415 {{0, 2, 1}},
1416 {{2, 0, 1}}}};
1417
1418 return table[face_orientation][line];
1419 }
1420 else if (*this == ReferenceCells::Pyramid)
1421 {
1422 if (face == 0) // The quadrilateral face
1423 {
1425 line,
1426 Utilities::get_bit(face_orientation, 0),
1427 Utilities::get_bit(face_orientation, 2),
1428 Utilities::get_bit(face_orientation, 1));
1429 }
1430 else // One of the triangular faces
1431 {
1432 static const ndarray<unsigned int, 6, 3> table = {{{{2, 1, 0}},
1433 {{0, 1, 2}},
1434 {{1, 0, 2}},
1435 {{1, 2, 0}},
1436 {{0, 2, 1}},
1437 {{2, 0, 1}}}};
1438
1439 return table[face_orientation][line];
1440 }
1441 }
1442 else if (*this == ReferenceCells::Wedge)
1443 {
1444 if (face > 1) // One of the quadrilateral faces
1445 {
1447 line,
1448 Utilities::get_bit(face_orientation, 0),
1449 Utilities::get_bit(face_orientation, 2),
1450 Utilities::get_bit(face_orientation, 1));
1451 }
1452 else // One of the triangular faces
1453 {
1454 static const ndarray<unsigned int, 6, 3> table = {{{{2, 1, 0}},
1455 {{0, 1, 2}},
1456 {{1, 0, 2}},
1457 {{1, 2, 0}},
1458 {{0, 2, 1}},
1459 {{2, 0, 1}}}};
1460
1461 return table[face_orientation][line];
1462 }
1463 }
1464 else if (*this == ReferenceCells::Hexahedron)
1465 {
1467 line,
1468 Utilities::get_bit(face_orientation, 0),
1469 Utilities::get_bit(face_orientation, 2),
1470 Utilities::get_bit(face_orientation, 1));
1471 }
1472
1473 Assert(false, ExcNotImplemented());
1475}
1476
1477
1478
1479namespace ReferenceCells
1480{
1481 template <int dim>
1482 inline constexpr const ReferenceCell &
1484 {
1485 switch (dim)
1486 {
1487 case 0:
1489 case 1:
1490 return ReferenceCells::Line;
1491 case 2:
1493 case 3:
1495 default:
1496 Assert(false, ExcNotImplemented());
1498 }
1499 }
1500
1501
1502
1503 template <int dim>
1504 inline constexpr const ReferenceCell &
1506 {
1507 switch (dim)
1508 {
1509 case 0:
1511 case 1:
1512 return ReferenceCells::Line;
1513 case 2:
1515 case 3:
1517 default:
1518 Assert(false, ExcNotImplemented());
1520 }
1521 }
1522} // namespace ReferenceCells
1523
1524
1525inline ReferenceCell
1526ReferenceCell::n_vertices_to_type(const int dim, const unsigned int n_vertices)
1527{
1528 AssertIndexRange(dim, 4);
1530
1531 const auto X = ReferenceCells::Invalid;
1532 static const ndarray<ReferenceCell, 4, 9> table = {
1533 {// dim 0
1534 {{X, ReferenceCells::Vertex, X, X, X, X, X, X, X}},
1535 // dim 1
1536 {{X, X, ReferenceCells::Line, X, X, X, X, X, X}},
1537 // dim 2
1538 {{X,
1539 X,
1540 X,
1543 X,
1544 X,
1545 X,
1546 X}},
1547 // dim 3
1548 {{X,
1549 X,
1550 X,
1551 X,
1555 X,
1558 ExcMessage("The combination of dim = " + std::to_string(dim) +
1559 " and n_vertices = " + std::to_string(n_vertices) +
1560 " does not correspond to a known reference cell type."));
1561 return table[dim][n_vertices];
1562}
1563
1564
1565
1566template <int dim>
1567inline double
1569 const unsigned int i) const
1570{
1572 if (*this == ReferenceCells::get_hypercube<dim>())
1574
1575 if (*this ==
1576 ReferenceCells::Triangle) // see also
1577 // BarycentricPolynomials<2>::compute_value
1578 {
1579 switch (i)
1580 {
1581 case 0:
1582 return 1.0 - xi[std::min(0, dim - 1)] - xi[std::min(1, dim - 1)];
1583 case 1:
1584 return xi[std::min(0, dim - 1)];
1585 case 2:
1586 return xi[std::min(1, dim - 1)];
1587 }
1588 }
1589
1590 if (*this ==
1591 ReferenceCells::Tetrahedron) // see also
1592 // BarycentricPolynomials<3>::compute_value
1593 {
1594 switch (i)
1595 {
1596 case 0:
1597 return 1.0 - xi[std::min(0, dim - 1)] - xi[std::min(1, dim - 1)] -
1598 xi[std::min(2, dim - 1)];
1599 case 1:
1600 return xi[std::min(0, dim - 1)];
1601 case 2:
1602 return xi[std::min(1, dim - 1)];
1603 case 3:
1604 return xi[std::min(2, dim - 1)];
1605 }
1606 }
1607
1608 if (*this ==
1609 ReferenceCells::Wedge) // see also
1610 // ScalarLagrangePolynomialWedge::compute_value
1611 {
1613 .d_linear_shape_function<2>(Point<2>(xi[std::min(0, dim - 1)],
1614 xi[std::min(1, dim - 1)]),
1615 i % 3) *
1617 .d_linear_shape_function<1>(Point<1>(xi[std::min(2, dim - 1)]),
1618 i / 3);
1619 }
1620
1621 if (*this ==
1622 ReferenceCells::Pyramid) // see also
1623 // ScalarLagrangePolynomialPyramid::compute_value
1624 {
1625 const double Q14 = 0.25;
1626 double ration;
1627
1628 const double r = xi[std::min(0, dim - 1)];
1629 const double s = xi[std::min(1, dim - 1)];
1630 const double t = xi[std::min(2, dim - 1)];
1631
1632 if (fabs(t - 1.0) > 1.0e-14)
1633 {
1634 ration = (r * s * t) / (1.0 - t);
1635 }
1636 else
1637 {
1638 ration = 0.0;
1639 }
1640
1641 if (i == 0)
1642 return Q14 * ((1.0 - r) * (1.0 - s) - t + ration);
1643 if (i == 1)
1644 return Q14 * ((1.0 + r) * (1.0 - s) - t - ration);
1645 if (i == 2)
1646 return Q14 * ((1.0 - r) * (1.0 + s) - t - ration);
1647 if (i == 3)
1648 return Q14 * ((1.0 + r) * (1.0 + s) - t + ration);
1649 else
1650 return t;
1651 }
1652
1653 Assert(false, ExcNotImplemented());
1654
1655 return 0.0;
1656}
1657
1658
1659
1660template <int dim>
1661inline Tensor<1, dim>
1663 const unsigned int i) const
1664{
1666 if (*this == ReferenceCells::get_hypercube<dim>())
1668
1669 if (*this ==
1670 ReferenceCells::Triangle) // see also
1671 // BarycentricPolynomials<2>::compute_grad
1672 {
1673 switch (i)
1674 {
1675 case 0:
1676 return Point<dim>(-1.0, -1.0);
1677 case 1:
1678 return Point<dim>(+1.0, +0.0);
1679 case 2:
1680 return Point<dim>(+0.0, +1.0);
1681 }
1682 }
1683
1684 Assert(false, ExcNotImplemented());
1685
1686 return Point<dim>(+0.0, +0.0, +0.0);
1687}
1688
1689
1690template <int dim>
1691inline Tensor<1, dim>
1693 const unsigned int i) const
1694{
1696 AssertIndexRange(i, dim - 1);
1697
1698 if (*this == ReferenceCells::get_hypercube<dim>())
1699 {
1702 }
1703 else if (*this == ReferenceCells::Triangle)
1704 {
1705 AssertIndexRange(face_no, 3);
1706 static const std::array<Tensor<1, dim>, 3> table = {
1707 {Point<dim>(1, 0),
1708 Point<dim>(-std::sqrt(0.5), +std::sqrt(0.5)),
1709 Point<dim>(0, -1)}};
1710
1711 return table[face_no];
1712 }
1713 else if (*this == ReferenceCells::Tetrahedron)
1714 {
1715 AssertIndexRange(face_no, 4);
1716 static const ndarray<Tensor<1, dim>, 4, 2> table = {
1717 {{{Point<dim>(0, 1, 0), Point<dim>(1, 0, 0)}},
1718 {{Point<dim>(1, 0, 0), Point<dim>(0, 0, 1)}},
1719 {{Point<dim>(0, 0, 1), Point<dim>(0, 1, 0)}},
1720 {{Point<dim>(-std::pow(1.0 / 3.0, 1.0 / 4.0),
1721 +std::pow(1.0 / 3.0, 1.0 / 4.0),
1722 0),
1723 Point<dim>(-std::pow(1.0 / 3.0, 1.0 / 4.0),
1724 0,
1725 +std::pow(1.0 / 3.0, 1.0 / 4.0))}}}};
1726
1727 return table[face_no][i];
1728 }
1729 else if (*this == ReferenceCells::Wedge)
1730 {
1731 AssertIndexRange(face_no, 5);
1732 static const ndarray<Tensor<1, dim>, 5, 2> table = {
1733 {{{Point<dim>(0, 1, 0), Point<dim>(1, 0, 0)}},
1734 {{Point<dim>(1, 0, 0), Point<dim>(0, 1, 0)}},
1735 {{Point<dim>(1, 0, 0), Point<dim>(0, 0, 1)}},
1736 {{Point<dim>(-1 / std::sqrt(2.0), +1 / std::sqrt(2.0), 0),
1737 Point<dim>(0, 0, 1)}},
1738 {{Point<dim>(0, 0, 1), Point<dim>(0, 1, 0)}}}};
1739
1740 return table[face_no][i];
1741 }
1742 else if (*this == ReferenceCells::Pyramid)
1743 {
1744 AssertIndexRange(face_no, 5);
1745 static const ndarray<Tensor<1, dim>, 5, 2> table = {
1746 {{{Point<dim>(0, 1, 0), Point<dim>(1, 0, 0)}},
1747 {{Point<dim>(+1.0 / sqrt(2.0), 0, +1.0 / sqrt(2.0)),
1748 Point<dim>(0, 1, 0)}},
1749 {{Point<dim>(+1.0 / sqrt(2.0), 0, -1.0 / sqrt(2.0)),
1750 Point<dim>(0, 1, 0)}},
1751 {{Point<dim>(1, 0, 0),
1752 Point<dim>(0, +1.0 / sqrt(2.0), +1.0 / sqrt(2.0))}},
1753 {{Point<dim>(1, 0, 0),
1754 Point<dim>(0, +1.0 / sqrt(2.0), -1.0 / sqrt(2.0))}}}};
1755
1756 return table[face_no][i];
1757 }
1758
1759 Assert(false, ExcNotImplemented());
1760
1761 return {};
1762}
1763
1764
1765
1766template <int dim>
1767inline Tensor<1, dim>
1768ReferenceCell::unit_normal_vectors(const unsigned int face_no) const
1769{
1770 AssertDimension(dim, this->get_dimension());
1771
1772 if (is_hyper_cube())
1773 {
1776 }
1777 else if (dim == 2)
1778 {
1780
1781 // Return the rotated vector
1782 return cross_product_2d(unit_tangential_vectors<dim>(face_no, 0));
1783 }
1784 else if (dim == 3)
1785 {
1786 return cross_product_3d(unit_tangential_vectors<dim>(face_no, 0),
1787 unit_tangential_vectors<dim>(face_no, 1));
1788 }
1789
1790 Assert(false, ExcNotImplemented());
1791
1792 return {};
1793}
1794
1795
1796
1797inline bool
1799 const unsigned int line,
1800 const unsigned char face_orientation_raw,
1801 const unsigned char line_orientation) const
1802{
1803 if (*this == ReferenceCells::Hexahedron)
1804 {
1805 static const bool bool_table[2][2][2][2] = {
1806 {{{true, false}, // lines 0/1, face_orientation=false,
1807 // face_flip=false, face_rotation=false and true
1808 {false, true}}, // lines 0/1, face_orientation=false,
1809 // face_flip=true, face_rotation=false and true
1810 {{true, true}, // lines 0/1, face_orientation=true,
1811 // face_flip=false, face_rotation=false and true
1812 {false, false}}}, // lines 0/1, face_orientation=true,
1813 // face_flip=true, face_rotation=false and true
1814
1815 {{{true, true}, // lines 2/3 ...
1816 {false, false}},
1817 {{true, false}, {false, true}}}};
1818
1819 const bool face_orientation = Utilities::get_bit(face_orientation_raw, 0);
1820 const bool face_flip = Utilities::get_bit(face_orientation_raw, 2);
1821 const bool face_rotation = Utilities::get_bit(face_orientation_raw, 1);
1822
1823 return (static_cast<bool>(line_orientation) ==
1824 bool_table[line / 2][face_orientation][face_flip][face_rotation]);
1825 }
1826 else
1827 // TODO: This might actually be wrong for some of the other
1828 // kinds of objects. We should check this
1829 return true;
1830}
1831
1832
1833
1834namespace internal
1835{
1836 template <typename T, std::size_t N>
1838 {
1839 public:
1843 NoPermutation(const ::ReferenceCell &entity_type,
1844 const std::array<T, N> & vertices_0,
1845 const std::array<T, N> & vertices_1)
1849 {}
1850
1854 virtual ~NoPermutation() noexcept override = default;
1855
1859 virtual void
1860 print_info(std::ostream &out) const override
1861 {
1862 out << "[";
1863
1864 const unsigned int n_vertices = entity_type.n_vertices();
1865
1866 for (unsigned int i = 0; i < n_vertices; ++i)
1867 {
1868 out << vertices_0[i];
1869 if (i + 1 != n_vertices)
1870 out << ",";
1871 }
1872
1873 out << "] is not a permutation of [";
1874
1875 for (unsigned int i = 0; i < n_vertices; ++i)
1876 {
1877 out << vertices_1[i];
1878 if (i + 1 != n_vertices)
1879 out << ",";
1880 }
1881
1882 out << "]." << std::endl;
1883 }
1884
1888 const ::ReferenceCell entity_type;
1889
1893 const std::array<T, N> vertices_0;
1894
1898 const std::array<T, N> vertices_1;
1899 };
1900} // namespace internal
1901
1902
1903
1904template <typename T, std::size_t N>
1905inline unsigned char
1906ReferenceCell::compute_orientation(const std::array<T, N> &vertices_0,
1907 const std::array<T, N> &vertices_1) const
1908{
1910 if (*this == ReferenceCells::Line)
1911 {
1912 const std::array<T, 2> i{{vertices_0[0], vertices_0[1]}};
1913 const std::array<T, 2> j{{vertices_1[0], vertices_1[1]}};
1914
1915 // line_orientation=true
1916 if (i == std::array<T, 2>{{j[0], j[1]}})
1917 return 1;
1918
1919 // line_orientation=false
1920 if (i == std::array<T, 2>{{j[1], j[0]}})
1921 return 0;
1922 }
1923 else if (*this == ReferenceCells::Triangle)
1924 {
1925 const std::array<T, 3> i{{vertices_0[0], vertices_0[1], vertices_0[2]}};
1926 const std::array<T, 3> j{{vertices_1[0], vertices_1[1], vertices_1[2]}};
1927
1928 // face_orientation=true, face_rotation=false, face_flip=false
1929 if (i == std::array<T, 3>{{j[0], j[1], j[2]}})
1930 return 1;
1931
1932 // face_orientation=true, face_rotation=true, face_flip=false
1933 if (i == std::array<T, 3>{{j[1], j[2], j[0]}})
1934 return 3;
1935
1936 // face_orientation=true, face_rotation=false, face_flip=true
1937 if (i == std::array<T, 3>{{j[2], j[0], j[1]}})
1938 return 5;
1939
1940 // face_orientation=false, face_rotation=false, face_flip=false
1941 if (i == std::array<T, 3>{{j[0], j[2], j[1]}})
1942 return 0;
1943
1944 // face_orientation=false, face_rotation=true, face_flip=false
1945 if (i == std::array<T, 3>{{j[2], j[1], j[0]}})
1946 return 2;
1947
1948 // face_orientation=false, face_rotation=false, face_flip=true
1949 if (i == std::array<T, 3>{{j[1], j[0], j[2]}})
1950 return 4;
1951 }
1952 else if (*this == ReferenceCells::Quadrilateral)
1953 {
1954 const std::array<T, 4> i{
1955 {vertices_0[0], vertices_0[1], vertices_0[2], vertices_0[3]}};
1956 const std::array<T, 4> j{
1957 {vertices_1[0], vertices_1[1], vertices_1[2], vertices_1[3]}};
1958
1959 // face_orientation=true, face_rotation=false, face_flip=false
1960 if (i == std::array<T, 4>{{j[0], j[1], j[2], j[3]}})
1961 return 1;
1962
1963 // face_orientation=true, face_rotation=true, face_flip=false
1964 if (i == std::array<T, 4>{{j[2], j[0], j[3], j[1]}})
1965 return 3;
1966
1967 // face_orientation=true, face_rotation=false, face_flip=true
1968 if (i == std::array<T, 4>{{j[3], j[2], j[1], j[0]}})
1969 return 5;
1970
1971 // face_orientation=true, face_rotation=true, face_flip=true
1972 if (i == std::array<T, 4>{{j[1], j[3], j[0], j[2]}})
1973 return 7;
1974
1975 // face_orientation=false, face_rotation=false, face_flip=false
1976 if (i == std::array<T, 4>{{j[0], j[2], j[1], j[3]}})
1977 return 0;
1978
1979 // face_orientation=false, face_rotation=true, face_flip=false
1980 if (i == std::array<T, 4>{{j[2], j[3], j[0], j[1]}})
1981 return 2;
1982
1983 // face_orientation=false, face_rotation=false, face_flip=true
1984 if (i == std::array<T, 4>{{j[3], j[1], j[2], j[0]}})
1985 return 4;
1986
1987 // face_orientation=false, face_rotation=true, face_flip=true
1988 if (i == std::array<T, 4>{{j[1], j[0], j[3], j[2]}})
1989 return 6;
1990 }
1991
1992 Assert(false, (internal::NoPermutation<T, N>(*this, vertices_0, vertices_1)));
1993
1994 return -1;
1995}
1996
1997
1998
1999template <typename T, std::size_t N>
2000inline std::array<T, N>
2002 const std::array<T, N> &vertices,
2003 const unsigned int orientation) const
2004{
2005 std::array<T, 4> temp;
2006
2007 if (*this == ReferenceCells::Line)
2008 {
2009 switch (orientation)
2010 {
2011 case 1:
2012 temp = {{vertices[0], vertices[1]}};
2013 break;
2014 case 0:
2015 temp = {{vertices[1], vertices[0]}};
2016 break;
2017 default:
2018 Assert(false, ExcNotImplemented());
2019 }
2020 }
2021 else if (*this == ReferenceCells::Triangle)
2022 {
2023 switch (orientation)
2024 {
2025 case 1:
2026 temp = {{vertices[0], vertices[1], vertices[2]}};
2027 break;
2028 case 3:
2029 temp = {{vertices[1], vertices[2], vertices[0]}};
2030 break;
2031 case 5:
2032 temp = {{vertices[2], vertices[0], vertices[1]}};
2033 break;
2034 case 0:
2035 temp = {{vertices[0], vertices[2], vertices[1]}};
2036 break;
2037 case 2:
2038 temp = {{vertices[2], vertices[1], vertices[0]}};
2039 break;
2040 case 4:
2041 temp = {{vertices[1], vertices[0], vertices[2]}};
2042 break;
2043 default:
2044 Assert(false, ExcNotImplemented());
2045 }
2046 }
2047 else if (*this == ReferenceCells::Quadrilateral)
2048 {
2049 switch (orientation)
2050 {
2051 case 1:
2052 temp = {{vertices[0], vertices[1], vertices[2], vertices[3]}};
2053 break;
2054 case 3:
2055 temp = {{vertices[2], vertices[0], vertices[3], vertices[1]}};
2056 break;
2057 case 5:
2058 temp = {{vertices[3], vertices[2], vertices[1], vertices[0]}};
2059 break;
2060 case 7:
2061 temp = {{vertices[1], vertices[3], vertices[0], vertices[2]}};
2062 break;
2063 case 0:
2064 temp = {{vertices[0], vertices[2], vertices[1], vertices[3]}};
2065 break;
2066 case 2:
2067 temp = {{vertices[2], vertices[3], vertices[0], vertices[1]}};
2068 break;
2069 case 4:
2070 temp = {{vertices[3], vertices[1], vertices[2], vertices[0]}};
2071 break;
2072 case 6:
2073 temp = {{vertices[1], vertices[0], vertices[3], vertices[2]}};
2074 break;
2075 default:
2076 Assert(false, ExcNotImplemented());
2077 }
2078 }
2079 else
2080 {
2082 }
2083
2084 std::array<T, N> temp_;
2085 std::copy_n(temp.begin(), N, temp_.begin());
2086
2087 return temp_;
2088}
2089
2090
2092
2093#endif
Abstract base class for mapping classes.
Definition: mapping.h:304
Definition: point.h:111
ArrayView< const unsigned int > faces_for_given_vertex(const unsigned int vertex_index) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
unsigned int n_vertices() const
bool is_hyper_cube() const
unsigned int vtk_linear_type() const
bool standard_vs_true_line_orientation(const unsigned int line, const unsigned char face_orientation, const unsigned char line_orientation) const
void serialize(Archive &archive, const unsigned int)
unsigned int child_cell_on_face(const unsigned int face_n, const unsigned int subface_n, const unsigned char face_orientation=1) const
unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const unsigned char face_orientation) const
unsigned char compute_orientation(const std::array< T, N > &vertices_0, const std::array< T, N > &vertices_1) const
std::array< unsigned int, 2 > standard_vertex_to_face_and_vertex_index(const unsigned int vertex) const
Quadrature< dim > get_gauss_type_quadrature(const unsigned n_points_1D) const
unsigned int standard_to_real_face_vertex(const unsigned int vertex, const unsigned int face, const unsigned char face_orientation) const
std::array< unsigned int, 2 > standard_line_to_face_and_line_index(const unsigned int line) const
double d_linear_shape_function(const Point< dim > &xi, const unsigned int i) const
Tensor< 1, dim > unit_tangential_vectors(const unsigned int face_no, const unsigned int i) const
std::array< T, N > permute_according_orientation(const std::array< T, N > &vertices, const unsigned int orientation) const
unsigned int n_faces() const
std::uint8_t kind
constexpr ReferenceCell()
unsigned int exodusii_vertex_to_deal_vertex(const unsigned int vertex_n) const
unsigned int unv_vertex_to_deal_vertex(const unsigned int vertex_n) const
std::unique_ptr< Mapping< dim, spacedim > > get_default_mapping(const unsigned int degree) const
unsigned int vtk_quadratic_type() const
unsigned int n_lines() const
unsigned int vtk_lagrange_type() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > line_indices() const
unsigned int get_dimension() const
ReferenceCell face_reference_cell(const unsigned int face_no) const
unsigned int exodusii_face_to_deal_face(const unsigned int face_n) const
static ReferenceCell n_vertices_to_type(const int dim, const unsigned int n_vertices)
unsigned int standard_to_real_face_line(const unsigned int line, const unsigned int face, const unsigned char face_orientation) const
const Mapping< dim, spacedim > & get_default_linear_mapping() const
std::string to_string() const
bool is_simplex() const
unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const unsigned char face_orientation) const
constexpr bool operator!=(const ReferenceCell &type) const
constexpr bool operator==(const ReferenceCell &type) const
const Quadrature< dim > & get_nodal_type_quadrature() const
Tensor< 1, dim > unit_normal_vectors(const unsigned int face_no) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices() const
Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i) const
const ::ReferenceCell entity_type
const std::array< T, N > vertices_0
const std::array< T, N > vertices_1
virtual void print_info(std::ostream &out) const override
virtual ~NoPermutation() noexcept override=default
NoPermutation(const ::ReferenceCell &entity_type, const std::array< T, N > &vertices_0, const std::array< T, N > &vertices_1)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_CONSTEXPR
Definition: config.h:172
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
Point< 3 > vertices[4]
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
std::string to_string(const T &t)
Definition: patterns.h:2329
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
Expression fabs(const Expression &x)
static const char U
static const char N
constexpr const ReferenceCell Tetrahedron
constexpr const ReferenceCell Quadrilateral
constexpr const ReferenceCell Wedge
constexpr const ReferenceCell Pyramid
constexpr const ReferenceCell Invalid
constexpr const ReferenceCell Triangle
constexpr const ReferenceCell Hexahedron
constexpr const ReferenceCell & get_hypercube()
constexpr const ReferenceCell Vertex
constexpr const ReferenceCell Line
constexpr const ReferenceCell & get_simplex()
bool get_bit(const unsigned char number, const unsigned int n)
Definition: utilities.h:1376
constexpr ::ReferenceCell make_reference_cell_from_int(const std::uint8_t kind)
static const unsigned int invalid_unsigned_int
Definition: types.h:196
boost::integer_range< IncrementableType > iota_view
Definition: iota_view.h:46
STL namespace.
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray
Definition: ndarray.h:108
static unsigned int standard_to_real_line_vertex(const unsigned int vertex, const bool line_orientation=true)
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim - 1 > &face_refinement_case=RefinementCase< dim - 1 >::isotropic_refinement)
static unsigned int standard_to_real_face_vertex(const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static std::array< unsigned int, 2 > standard_hex_line_to_quad_line_index(const unsigned int line)
static unsigned int standard_to_real_face_line(const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static std::array< unsigned int, 2 > standard_hex_vertex_to_quad_vertex_index(const unsigned int vertex)
static std::array< unsigned int, 2 > standard_quad_vertex_to_line_vertex_index(const unsigned int vertex)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
Definition: tensor.h:2539
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
Definition: tensor.h:2564