Reference documentation for deal.II version 9.6.0
\(\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\}}\)
Loading...
Searching...
No Matches
reference_cell.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2020 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_tria_reference_cell_h
16#define dealii_tria_reference_cell_h
17
18#include <deal.II/base/config.h>
19
23#include <deal.II/base/point.h>
24#include <deal.II/base/tensor.h>
26
28
29#include <boost/container/small_vector.hpp>
30
31#include <iosfwd>
32#include <string>
33
35
36// Forward declarations
37#ifndef DOXYGEN
38template <int dim, int spacedim>
39class Mapping;
40
41template <int dim>
42class Quadrature;
43
44class ReferenceCell;
45#endif
46
47
48namespace internal
49{
62 constexpr ReferenceCell
63 make_reference_cell_from_int(const std::uint8_t kind);
64} // namespace internal
65
66
67
74enum class IsotropicRefinementChoice : std::uint8_t
75{
83 cut_tet_68 = 1,
87 cut_tet_57 = 2,
91 cut_tet_49 = 3,
92};
93
94
95
127{
128public:
136 static ReferenceCell
137 n_vertices_to_type(const int dim, const unsigned int n_vertices);
138
147 constexpr ReferenceCell();
148
159 bool
160 is_hyper_cube() const;
161
165 bool
166 is_simplex() const;
167
172 unsigned int
173 get_dimension() const;
174
188 template <int dim>
189 double
190 d_linear_shape_function(const Point<dim> &xi, const unsigned int i) const;
191
196 template <int dim>
199 const unsigned int i) const;
200
210 template <int dim, int spacedim = dim>
211 std::unique_ptr<Mapping<dim, spacedim>>
212 get_default_mapping(const unsigned int degree) const;
213
225 template <int dim, int spacedim = dim>
228
237 template <int dim>
239 get_gauss_type_quadrature(const unsigned n_points_1d) const;
240
250 template <int dim>
253
267 template <int dim>
268 const Quadrature<dim> &
270
285 unsigned int
286 n_vertices() const;
287
293 vertex_indices() const;
294
302 template <int dim>
304 vertex(const unsigned int v) const;
305
311 unsigned int
312 n_lines() const;
313
319 line_indices() const;
320
326 unsigned int
327 n_faces() const;
328
334 face_indices() const;
335
339 unsigned int
341
346 get_isotropic_refinement_choice(const unsigned int ref_choice) const;
347
359 unsigned int
360 n_isotropic_children() const;
361
366 template <int dim>
367 unsigned int
368 n_children(const RefinementCase<dim> ref_case =
370
377
391 face_reference_cell(const unsigned int face_no) const;
392
406 static constexpr unsigned char
408
419 static constexpr unsigned char
421
463 unsigned int
464 child_cell_on_face(const unsigned int face,
465 const unsigned int subface,
466 const unsigned char face_orientation =
468
478 std::array<unsigned int, 2>
479 standard_vertex_to_face_and_vertex_index(const unsigned int vertex) const;
480
490 std::array<unsigned int, 2>
491 standard_line_to_face_and_line_index(const unsigned int line) const;
492
501 unsigned int
502 line_to_cell_vertices(const unsigned int line,
503 const unsigned int vertex) const;
504
508 unsigned int
509 face_to_cell_lines(const unsigned int face,
510 const unsigned int line,
511 const unsigned char face_orientation) const;
512
516 unsigned int
517 face_to_cell_vertices(const unsigned int face,
518 const unsigned int vertex,
519 const unsigned char face_orientation) const;
520
540 template <int dim>
542 face_vertex_location(const unsigned int face,
543 const unsigned int vertex) const;
544
548 unsigned int
549 standard_to_real_face_vertex(const unsigned int vertex,
550 const unsigned int face,
551 const unsigned char face_orientation) const;
552
556 unsigned int
557 standard_to_real_face_line(const unsigned int line,
558 const unsigned int face,
559 const unsigned char face_orientation) const;
560
571 bool
572 standard_vs_true_line_orientation(const unsigned int line,
573 const unsigned int face,
574 const unsigned char face_orientation,
575 const bool line_orientation) const;
576
600 double
601 volume() const;
602
615 template <int dim>
617 barycenter() const;
618
638 template <int dim>
639 bool
640 contains_point(const Point<dim> &p, const double tolerance = 0) const;
641
646 template <int dim>
648 closest_point(const Point<dim> &p) const;
649
657 template <int dim>
659 unit_tangential_vectors(const unsigned int face_no,
660 const unsigned int i) const;
661
665 template <int dim>
667 unit_normal_vectors(const unsigned int face_no) const;
668
674 unsigned int
675 n_face_orientations(const unsigned int face_no) const;
676
693 template <typename T, std::size_t N>
694 DEAL_II_DEPRECATED unsigned char
695 compute_orientation(const std::array<T, N> &vertices_0,
696 const std::array<T, N> &vertices_1) const;
697
723 template <typename T>
724 unsigned char
726 const ArrayView<const T> &vertices_1) const;
727
728
742 template <typename T, std::size_t N>
743 DEAL_II_DEPRECATED std::array<T, N>
744 permute_according_orientation(const std::array<T, N> &vertices,
745 const unsigned int orientation) const;
746
758 template <typename T>
759 boost::container::small_vector<T, 8>
761 const unsigned char orientation) const;
762
768 unsigned char
769 get_inverse_combined_orientation(const unsigned char orientation) const;
770
775 faces_for_given_vertex(const unsigned int vertex_index) const;
776
781 constexpr ::ndarray<unsigned int, 12, 4>
782 new_isotropic_child_face_lines(const unsigned int refinement_choice) const;
783
788 constexpr ::ndarray<unsigned int, 12, 4, 2>
790 const unsigned int refinement_choice) const;
791
796 constexpr ::ndarray<unsigned int, 8, 6>
797 new_isotropic_child_cell_faces(const unsigned int refinement_choice) const;
798
803 constexpr ::ndarray<unsigned int, 8, 4>
804 new_isotropic_child_cell_vertices(const unsigned int refinement_choice) const;
805
818 unsigned int
819 exodusii_vertex_to_deal_vertex(const unsigned int vertex_n) const;
820
824 unsigned int
825 exodusii_face_to_deal_face(const unsigned int face_n) const;
826
830 unsigned int
831 unv_vertex_to_deal_vertex(const unsigned int vertex_n) const;
832
836 unsigned int
837 vtk_linear_type() const;
838
843 unsigned int
844 vtk_quadratic_type() const;
845
850 unsigned int
851 vtk_lagrange_type() const;
852
865 template <int dim>
866 unsigned int
868 const std::array<unsigned, dim> &node_indices,
869 const std::array<unsigned, dim> &nodes_per_direction,
870 const bool legacy_format) const;
871
875 unsigned int
876 vtk_vertex_to_deal_vertex(const unsigned int vertex_index) const;
877
881 unsigned int
882 gmsh_element_type() const;
883
897 std::string
898 to_string() const;
899
903 constexpr operator std::uint8_t() const;
904
908 constexpr bool
909 operator==(const ReferenceCell &type) const;
910
914 constexpr bool
915 operator!=(const ReferenceCell &type) const;
916
922 template <class Archive>
923 void
924 serialize(Archive &archive, const unsigned int /*version*/);
925
929 static constexpr std::size_t
931
936private:
940 std::uint8_t kind;
941
948 constexpr ReferenceCell(const std::uint8_t kind);
949
956 {{{0}}}};
957
962 {{{1, 0}}, {{0, 1}}}};
963
968 {{{0, 2, 1}},
969 {{0, 1, 2}},
970 {{2, 1, 0}},
971 {{2, 0, 1}},
972 {{1, 0, 2}},
973 {{1, 2, 0}}}};
974
978 static constexpr ndarray<unsigned int, 8, 4>
980 {{0, 1, 2, 3}},
981 {{2, 3, 0, 1}},
982 {{2, 0, 3, 1}},
983 {{3, 1, 2, 0}},
984 {{3, 2, 1, 0}},
985 {{1, 0, 3, 2}},
986 {{1, 3, 0, 2}}}};
987
992 friend constexpr ReferenceCell
994
995 friend std::ostream &
996 operator<<(std::ostream &out, const ReferenceCell &reference_cell);
997
998 friend std::istream &
999 operator>>(std::istream &in, ReferenceCell &reference_cell);
1000};
1001
1002
1010std::ostream &
1011operator<<(std::ostream &out, const ReferenceCell &reference_cell);
1012
1020std::istream &
1021operator>>(std::istream &in, ReferenceCell &reference_cell);
1022
1023
1024inline constexpr ReferenceCell::ReferenceCell(const std::uint8_t kind)
1025 : kind(kind)
1026{}
1027
1028
1029
1030inline constexpr ReferenceCell::operator std::uint8_t() const
1031{
1032 return kind;
1033}
1034
1035
1036
1037inline constexpr bool
1039{
1040 return kind == type.kind;
1041}
1042
1043
1044
1045inline constexpr bool
1047{
1048 return kind != type.kind;
1049}
1050
1051
1052
1053namespace internal
1054{
1055 inline constexpr ReferenceCell
1056 make_reference_cell_from_int(const std::uint8_t kind)
1057 {
1058#ifndef DEAL_II_CXX14_CONSTEXPR_BUG
1059 // Make sure these are the only indices from which objects can be
1060 // created.
1061 Assert((kind == std::numeric_limits<std::uint8_t>::max()) || (kind < 8),
1063#endif
1064
1065 // Call the private constructor, which we can from here because this
1066 // function is a 'friend'.
1067 return {kind};
1068 }
1069} // namespace internal
1070
1071
1072
1081{
1082 constexpr const ReferenceCell Vertex =
1084 constexpr const ReferenceCell Line =
1086 constexpr const ReferenceCell Triangle =
1090 constexpr const ReferenceCell Tetrahedron =
1092 constexpr const ReferenceCell Pyramid =
1094 constexpr const ReferenceCell Wedge =
1096 constexpr const ReferenceCell Hexahedron =
1098 constexpr const ReferenceCell Invalid =
1100 std::numeric_limits<std::uint8_t>::max());
1101
1107 template <int dim>
1108 constexpr const ReferenceCell &
1109 get_simplex();
1110
1116 template <int dim>
1117 constexpr const ReferenceCell &
1118 get_hypercube();
1119} // namespace ReferenceCells
1120
1121
1122
1124 : ReferenceCell(ReferenceCells::Invalid)
1125{}
1126
1127
1128
1129template <class Archive>
1130inline void
1131ReferenceCell::serialize(Archive &archive, const unsigned int /*version*/)
1132{
1133 archive &kind;
1134}
1135
1136
1137
1138inline constexpr std::size_t
1143
1144
1145
1147ReferenceCell::faces_for_given_vertex(const unsigned int vertex) const
1148{
1150 switch (this->kind)
1151 {
1153 return {&GeometryInfo<2>::vertex_to_face[vertex][0], 1};
1155 return {&GeometryInfo<2>::vertex_to_face[vertex][0], 2};
1157 {
1158 static constexpr ndarray<unsigned int, 3, 2> table = {
1159 {{{0, 2}}, {{0, 1}}, {{1, 2}}}};
1160 return table[vertex];
1161 }
1163 {
1164 static constexpr ndarray<unsigned int, 4, 3> table = {
1165 {{{0, 1, 2}}, {{0, 1, 3}}, {{0, 2, 3}}, {{1, 2, 3}}}};
1166
1167 return table[vertex];
1168 }
1170 {
1171 static constexpr unsigned int X = numbers::invalid_unsigned_int;
1172 static constexpr ndarray<unsigned int, 5, 4> table = {
1173 {{{0, 1, 3, X}},
1174 {{0, 2, 3, X}},
1175 {{0, 1, 4, X}},
1176 {{0, 2, 4, X}},
1177 {{1, 2, 3, 4}}}};
1178
1179 return {&table[vertex][0], vertex == 4 ? 4u : 3u};
1180 }
1182 {
1184 static constexpr ndarray<unsigned int, 6, 3> table = {{{{0, 2, 4}},
1185 {{0, 2, 3}},
1186 {{0, 3, 4}},
1187 {{1, 2, 4}},
1188 {{1, 2, 3}},
1189 {{1, 3, 4}}}};
1190
1191 return table[vertex];
1192 }
1194 return {&GeometryInfo<3>::vertex_to_face[vertex][0], 3};
1195 default:
1197 }
1198
1199 return {};
1200}
1201
1202
1203constexpr ::ndarray<unsigned int, 12, 4>
1205 const unsigned int refienement_choice) const
1206{
1207 // AssertIndexRange(refienement_choice, n_isotropic_refinement_choices());
1208 const unsigned int X = numbers::invalid_unsigned_int;
1209 switch (this->kind)
1210 {
1212 {
1213 constexpr ::ndarray<unsigned int, 3, 12, 4> new_quad_lines_tet =
1214 {{// new line is (6,8)
1215 {{{{2, 3, 8, X}},
1216 {{0, 9, 5, X}},
1217 {{1, 6, 11, X}},
1218 {{4, 10, 7, X}},
1219 {{2, 12, 5, X}},
1220 {{1, 9, 12, X}},
1221 {{4, 8, 12, X}},
1222 {{6, 12, 10, X}},
1223 {{X, X, X, X}},
1224 {{X, X, X, X}},
1225 {{X, X, X, X}},
1226 {{X, X, X, X}}}},
1227 // new line is (5,7)
1228 {{{{2, 3, 8, X}},
1229 {{0, 9, 5, X}},
1230 {{1, 6, 11, X}},
1231 {{4, 10, 7, X}},
1232 {{0, 3, 12, X}},
1233 {{1, 12, 8, X}},
1234 {{4, 12, 9, X}},
1235 {{7, 11, 12, X}},
1236 {{X, X, X, X}},
1237 {{X, X, X, X}},
1238 {{X, X, X, X}},
1239 {{X, X, X, X}}}},
1240 // new line is (4,9)
1241 {{{{2, 3, 8, X}},
1242 {{0, 9, 5, X}},
1243 {{1, 6, 11, X}},
1244 {{4, 10, 7, X}},
1245 {{0, 12, 11, X}},
1246 {{2, 6, 12, X}},
1247 {{3, 12, 7, X}},
1248 {{5, 10, 12, X}},
1249 {{X, X, X, X}},
1250 {{X, X, X, X}},
1251 {{X, X, X, X}},
1252 {{X, X, X, X}}}}}};
1253 return new_quad_lines_tet[refienement_choice];
1254 }
1255
1257 {
1258 constexpr ::ndarray<unsigned int, 12, 4> new_quad_lines_hex = {
1259 {{{10, 28, 16, 24}},
1260 {{28, 14, 17, 25}},
1261 {{11, 29, 24, 20}},
1262 {{29, 15, 25, 21}},
1263 {{18, 26, 0, 28}},
1264 {{26, 22, 1, 29}},
1265 {{19, 27, 28, 4}},
1266 {{27, 23, 29, 5}},
1267 {{2, 24, 8, 26}},
1268 {{24, 6, 9, 27}},
1269 {{3, 25, 26, 12}},
1270 {{25, 7, 27, 13}}}};
1271 return new_quad_lines_hex;
1272 }
1273 default:
1275 }
1276
1277 return {};
1278}
1279
1280
1281
1282constexpr ::ndarray<unsigned int, 12, 4, 2>
1284 const unsigned int refienement_choice) const
1285{
1286 // AssertIndexRange(refienement_choice, n_isotropic_refinement_choices());
1287 const unsigned int X = numbers::invalid_unsigned_int;
1288 switch (this->kind)
1289 {
1291 {
1292 constexpr ::ndarray<unsigned int, 3, 12, 4, 2> table_tet = {
1293 {// new line is (6, 8)
1294 {{{{{{6, 4}}, {{4, 7}}, {{7, 6}}, {{X, X}}}},
1295 {{{{4, 5}}, {{5, 8}}, {{8, 4}}, {{X, X}}}},
1296 {{{{5, 6}}, {{6, 9}}, {{9, 5}}, {{X, X}}}},
1297 {{{{7, 8}}, {{8, 9}}, {{9, 7}}, {{X, X}}}},
1298 {{{{4, 6}}, {{6, 8}}, {{8, 4}}, {{X, X}}}},
1299 {{{{6, 5}}, {{5, 8}}, {{8, 6}}, {{X, X}}}},
1300 {{{{8, 7}}, {{7, 6}}, {{6, 8}}, {{X, X}}}},
1301 {{{{9, 6}}, {{6, 8}}, {{8, 9}}, {{X, X}}}},
1302 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1303 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1304 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1305 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}}}},
1306 // new line is (5, 7)
1307 {{{{{{6, 4}}, {{4, 7}}, {{7, 6}}, {{X, X}}}},
1308 {{{{4, 5}}, {{5, 8}}, {{8, 4}}, {{X, X}}}},
1309 {{{{5, 6}}, {{6, 9}}, {{9, 5}}, {{X, X}}}},
1310 {{{{7, 8}}, {{8, 9}}, {{9, 7}}, {{X, X}}}},
1311 {{{{5, 4}}, {{4, 7}}, {{7, 5}}, {{X, X}}}},
1312 {{{{6, 5}}, {{5, 7}}, {{7, 6}}, {{X, X}}}},
1313 {{{{8, 7}}, {{7, 5}}, {{5, 8}}, {{X, X}}}},
1314 {{{{7, 9}}, {{9, 5}}, {{5, 7}}, {{X, X}}}},
1315 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1316 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1317 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1318 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}}}},
1319 // new line is (4, 9)
1320 {{{{{{6, 4}}, {{4, 7}}, {{7, 6}}, {{X, X}}}},
1321 {{{{4, 5}}, {{5, 8}}, {{8, 4}}, {{X, X}}}},
1322 {{{{5, 6}}, {{6, 9}}, {{9, 5}}, {{X, X}}}},
1323 {{{{7, 8}}, {{8, 9}}, {{9, 7}}, {{X, X}}}},
1324 {{{{5, 4}}, {{4, 9}}, {{9, 5}}, {{X, X}}}},
1325 {{{{4, 6}}, {{6, 9}}, {{9, 4}}, {{X, X}}}},
1326 {{{{7, 4}}, {{4, 9}}, {{9, 7}}, {{X, X}}}},
1327 {{{{4, 8}}, {{8, 9}}, {{9, 4}}, {{X, X}}}},
1328 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1329 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1330 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}},
1331 {{{{X, X}}, {{X, X}}, {{X, X}}, {{X, X}}}}}}}};
1332 return table_tet[refienement_choice];
1333 }
1334
1336 {
1337 constexpr ::ndarray<unsigned int, 12, 4, 2> table_hex = {
1338 {{{{{10, 22}}, {{24, 26}}, {{10, 24}}, {{22, 26}}}},
1339 {{{{24, 26}}, {{11, 23}}, {{24, 11}}, {{26, 23}}}},
1340 {{{{22, 14}}, {{26, 25}}, {{22, 26}}, {{14, 25}}}},
1341 {{{{26, 25}}, {{23, 15}}, {{26, 23}}, {{25, 15}}}},
1342 {{{{8, 24}}, {{20, 26}}, {{8, 20}}, {{24, 26}}}},
1343 {{{{20, 26}}, {{12, 25}}, {{20, 12}}, {{26, 25}}}},
1344 {{{{24, 9}}, {{26, 21}}, {{24, 26}}, {{9, 21}}}},
1345 {{{{26, 21}}, {{25, 13}}, {{26, 25}}, {{21, 13}}}},
1346 {{{{16, 20}}, {{22, 26}}, {{16, 22}}, {{20, 26}}}},
1347 {{{{22, 26}}, {{17, 21}}, {{22, 17}}, {{26, 21}}}},
1348 {{{{20, 18}}, {{26, 23}}, {{20, 26}}, {{18, 23}}}},
1349 {{{{26, 23}}, {{21, 19}}, {{26, 21}}, {{23, 19}}}}}};
1350 return table_hex;
1351 }
1352 default:
1354 }
1355
1356 return {};
1357}
1358
1359
1360
1361constexpr ::ndarray<unsigned int, 8, 6>
1363 const unsigned int refienement_choice) const
1364{
1365 // AssertIndexRange(refienement_choice, n_isotropic_refinement_choices());
1366 const unsigned int X = numbers::invalid_unsigned_int;
1367 switch (this->kind)
1368 {
1370 {
1371 constexpr ::ndarray<unsigned int, 3, 8, 6> cell_quads_tet = {
1372 {// new line is (6, 8)
1373 {{{{8, 13, 16, 0, X, X}},
1374 {{9, 12, 1, 21, X, X}},
1375 {{10, 2, 17, 20, X, X}},
1376 {{3, 14, 18, 22, X, X}},
1377 {{11, 1, 4, 5, X, X}},
1378 {{15, 0, 4, 6, X, X}},
1379 {{19, 7, 6, 3, X, X}},
1380 {{23, 5, 2, 7, X, X}}}},
1381 // new line is (5, 7)
1382 {{
1383 {{8, 13, 16, 0, X, X}},
1384 {{9, 12, 1, 21, X, X}},
1385 {{10, 2, 17, 20, X, X}},
1386 {{3, 14, 18, 22, X, X}},
1387 {{11, 4, 0, 5, X, X}},
1388 {{15, 4, 1, 6, X, X}},
1389 {{19, 2, 5, 7, X, X}},
1390 {{23, 6, 7, 3, X, X}},
1391 }},
1392 // new line is (4, 9)
1393 {{
1394 {{8, 13, 16, 0, X, X}},
1395 {{9, 12, 1, 21, X, X}},
1396 {{10, 2, 17, 20, X, X}},
1397 {{3, 14, 18, 22, X, X}},
1398 {{11, 4, 5, 2, X, X}},
1399 {{15, 6, 7, 3, X, X}},
1400 {{19, 5, 0, 6, X, X}},
1401 {{23, 1, 4, 7, X, X}},
1402 }}}};
1403 return cell_quads_tet[refienement_choice];
1404 }
1405
1407 {
1408 constexpr ::ndarray<unsigned int, 8, 6> cell_quads_hex = {{
1409 {{12, 0, 20, 4, 28, 8}}, // bottom children
1410 {{0, 16, 22, 6, 29, 9}}, //
1411 {{13, 1, 4, 24, 30, 10}}, //
1412 {{1, 17, 6, 26, 31, 11}}, //
1413 {{14, 2, 21, 5, 8, 32}}, // top children
1414 {{2, 18, 23, 7, 9, 33}}, //
1415 {{15, 3, 5, 25, 10, 34}}, //
1416 {{3, 19, 7, 27, 11, 35}} //
1417 }};
1418 return cell_quads_hex;
1419 }
1420 default:
1422 }
1423
1424 return {};
1425}
1426
1427
1428
1429constexpr ::ndarray<unsigned int, 8, 4>
1431 const unsigned int refienement_choice) const
1432{
1433 // AssertIndexRange(refienement_choice, n_isotropic_refinement_choices());
1434
1435 switch (this->kind)
1436 {
1438 {
1439 constexpr ::ndarray<unsigned int, 3, 8, 4> cell_vertices_tet = {
1440 {// new line is (6,8)
1441 {{{{0, 4, 6, 7}},
1442 {{4, 1, 5, 8}},
1443 {{6, 5, 2, 9}},
1444 {{7, 8, 9, 3}},
1445 {{4, 5, 6, 8}},
1446 {{4, 7, 8, 6}},
1447 {{6, 9, 7, 8}},
1448 {{5, 8, 9, 6}}}},
1449 // new line is (5,7)
1450 {{{{0, 4, 6, 7}},
1451 {{4, 1, 5, 8}},
1452 {{6, 5, 2, 9}},
1453 {{7, 8, 9, 3}},
1454 {{4, 5, 6, 7}},
1455 {{4, 7, 8, 5}},
1456 {{6, 9, 7, 5}},
1457 {{5, 8, 9, 7}}}},
1458 // new line is (4,9)
1459 {{{{0, 4, 6, 7}},
1460 {{4, 1, 5, 8}},
1461 {{6, 5, 2, 9}},
1462 {{7, 8, 9, 3}},
1463 {{4, 5, 6, 9}},
1464 {{4, 7, 8, 9}},
1465 {{6, 9, 7, 4}},
1466 {{5, 8, 9, 4}}}}}};
1467 return cell_vertices_tet[refienement_choice];
1468 }
1469 default:
1471 }
1472
1473 return {};
1474}
1475
1476
1477
1478inline bool
1480{
1481 return (*this == ReferenceCells::Vertex || *this == ReferenceCells::Line ||
1484}
1485
1486
1487
1488inline bool
1490{
1491 return (*this == ReferenceCells::Vertex || *this == ReferenceCells::Line ||
1492 *this == ReferenceCells::Triangle ||
1494}
1495
1496
1497
1498inline unsigned int
1500{
1501 switch (this->kind)
1502 {
1504 return 0;
1506 return 1;
1509 return 2;
1514 return 3;
1515 default:
1517 }
1518
1520}
1521
1522
1523
1524template <int dim>
1527{
1528 return Quadrature<dim>(std::vector<Point<dim>>({barycenter<dim>()}),
1529 std::vector<double>({volume()}));
1530}
1531
1532
1533
1534inline unsigned int
1536{
1537 switch (this->kind)
1538 {
1540 return 1;
1542 return 2;
1544 return 3;
1546 return 4;
1548 return 4;
1550 return 5;
1552 return 6;
1554 return 8;
1555 default:
1557 }
1558
1560}
1561
1562
1563
1564inline unsigned int
1566{
1567 switch (this->kind)
1568 {
1570 return 0;
1572 return 1;
1574 return 3;
1576 return 4;
1578 return 6;
1580 return 8;
1582 return 9;
1584 return 12;
1585 default:
1587 }
1588
1590}
1591
1592
1593
1594template <int dim>
1596ReferenceCell::vertex(const unsigned int v) const
1597{
1600
1601 switch (dim)
1602 {
1603 case 0:
1604 {
1605 if (*this == ReferenceCells::Vertex)
1606 return Point<dim>();
1607 break;
1608 }
1609 case 1:
1610 {
1611 static const Point<dim> vertices[2] = {
1612 Point<dim>(), // the origin
1613 Point<dim>::unit_vector(0) // unit point along x-axis
1614 };
1615 if (*this == ReferenceCells::Line)
1616 return vertices[v];
1617 break;
1618 }
1619 case 2:
1620 {
1621 switch (this->kind)
1622 {
1624 {
1625 static const Point<dim> vertices[3] = {
1626 Point<dim>(), // the origin
1627 Point<dim>::unit_vector(0), // unit point along x-axis
1628 Point<dim>::unit_vector(1) // unit point along y-axis
1629 };
1630 return vertices[v];
1631 }
1633 {
1634 static const Point<dim> vertices[4] = {
1635 // First the two points on the x-axis
1636 Point<dim>(),
1638 // Then these two points shifted in the y-direction
1641 return vertices[v];
1642 }
1643 }
1644 break;
1645 }
1646 case 3:
1647 {
1648 switch (this->kind)
1649 {
1651 {
1652 static const Point<dim> vertices[4] = {
1653 Point<dim>(), // the origin
1654 Point<dim>::unit_vector(0), // unit point along x-axis
1655 Point<dim>::unit_vector(1), // unit point along y-axis
1656 Point<dim>::unit_vector(2) // unit point along z-axis
1657 };
1658 return vertices[v];
1659 }
1661 {
1662 static const Point<dim> vertices[5] = {
1663 Point<dim>{-1.0, -1.0, 0.0},
1664 Point<dim>{+1.0, -1.0, 0.0},
1665 Point<dim>{-1.0, +1.0, 0.0},
1666 Point<dim>{+1.0, +1.0, 0.0},
1667 Point<dim>{+0.0, +0.0, 1.0}};
1668 return vertices[v];
1669 }
1671 {
1672 static const Point<dim> vertices[6] = {
1673 // First the three points on the triangular base of the
1674 // wedge:
1675 Point<dim>(),
1678 // And now everything shifted in the z-direction again
1682 return vertices[v];
1683 }
1685 {
1686 static const Point<dim> vertices[8] = {
1687 // First the two points on the x-axis
1688 Point<dim>(),
1690 // Then these two points shifted in the y-direction
1693 // And now all four points shifted in the z-direction
1700 return vertices[v];
1701 }
1702 }
1703 break;
1704 }
1705 default:
1707 }
1708
1710 return Point<dim>();
1711}
1712
1713
1714inline unsigned int
1716{
1717 switch (this->kind)
1718 {
1720 return 0;
1722 return 2;
1724 return 3;
1726 return 4;
1728 return 4;
1730 return 5;
1732 return 5;
1734 return 6;
1735 default:
1737 }
1738
1740}
1741
1742
1743
1746{
1747 return {0U, n_faces()};
1748}
1749
1750
1751
1752inline unsigned int
1754{
1755 switch (this->kind)
1756 {
1758 return 1;
1760 return 1;
1762 return 1;
1764 return 1;
1766 return 3;
1768 return 1;
1770 return 1;
1772 return 1;
1773 default:
1775 }
1776
1778}
1779
1780
1781
1784 const unsigned int ref_choice) const
1785{
1787 switch (this->kind)
1788 {
1790 {
1792 isotropic_ref_choices = {{IsotropicRefinementChoice::cut_tet_68,
1795 return isotropic_ref_choices[ref_choice];
1796 }
1797 default:
1798 {
1800 }
1801 }
1802
1804}
1805
1806
1807
1808inline unsigned int
1810{
1811 switch (this->kind)
1812 {
1814 return 0;
1816 return 2;
1818 return 4;
1820 return 4;
1822 return 8;
1824 // We haven't yet decided how to refine pyramids. Update this when we
1825 // have
1826 return 0;
1828 return 8;
1830 return 8;
1831 default:
1833 }
1834
1836}
1837
1838
1839
1840template <int dim>
1841unsigned int
1843{
1844 // Use GeometryInfo here to keep it the single source of truth
1845 if (this->is_hyper_cube())
1846 return GeometryInfo<dim>::n_children(ref_case);
1847 else
1848 return this->n_isotropic_children();
1849}
1850
1851
1852
1855{
1856 return {0U, n_isotropic_children()};
1857}
1858
1859
1860
1863{
1864 return {0U, n_vertices()};
1865}
1866
1867
1868
1871{
1872 return {0U, n_lines()};
1873}
1874
1875
1876
1877inline ReferenceCell
1878ReferenceCell::face_reference_cell(const unsigned int face_no) const
1879{
1880 AssertIndexRange(face_no, n_faces());
1881
1882 switch (this->kind)
1883 {
1890 return ReferenceCells::Line;
1894 if (face_no == 0)
1896 else
1899 if (face_no > 1)
1901 else
1905 default:
1907 }
1908
1910}
1911
1912
1913
1914inline constexpr unsigned char
1916{
1917 // Our convention is that 'orientation' has a default value of true and
1918 // occupies the least-significant bit while rotate and flip have default
1919 // values of 'false' and occupy the second and third bits.
1920 return 0b001;
1921}
1922
1923
1924
1925inline constexpr unsigned char
1927{
1928 // For a reversed line 'orientation' is false and neither flip nor rotate are
1929 // defined.
1930 return 0b000;
1931}
1932
1933
1934
1935inline unsigned int
1937 const unsigned int face,
1938 const unsigned int subface,
1939 const unsigned char combined_face_orientation) const
1940{
1941 AssertIndexRange(face, n_faces());
1943
1944 switch (this->kind)
1945 {
1948 {
1950 break;
1951 }
1953 {
1954 static constexpr ndarray<unsigned int, 3, 2> subcells = {
1955 {{{0, 1}}, {{1, 2}}, {{2, 0}}}};
1956
1957 return subcells[face][subface];
1958 }
1960 {
1961 const auto [face_orientation, face_rotation, face_flip] =
1962 internal::split_face_orientation(combined_face_orientation);
1963
1966 face,
1967 subface,
1968 face_orientation,
1969 face_flip,
1970 face_rotation);
1971 }
1975 {
1977 break;
1978 }
1980 {
1981 const auto [face_orientation, face_rotation, face_flip] =
1982 internal::split_face_orientation(combined_face_orientation);
1983
1986 face,
1987 subface,
1988 face_orientation,
1989 face_flip,
1990 face_rotation);
1991 }
1992 default:
1994 }
1995
1997}
1998
1999
2000
2001inline std::array<unsigned int, 2>
2003 const unsigned int vertex) const
2004{
2006 // Work around a GCC warning at higher optimization levels by making all of
2007 // these tables the same size
2008 constexpr unsigned int X = numbers::invalid_unsigned_int;
2009
2010 switch (this->kind)
2011 {
2015 break;
2017 {
2018 static constexpr ndarray<unsigned int, 6, 2> table = {
2019 {{{0, 0}}, {{0, 1}}, {{1, 1}}, {{X, X}}, {{X, X}}, {{X, X}}}};
2020
2021 return table[vertex];
2022 }
2024 {
2026 vertex);
2027 }
2029 {
2030 static constexpr ndarray<unsigned int, 6, 2> table = {
2031 {{{0, 0}}, {{0, 1}}, {{0, 2}}, {{1, 2}}, {{X, X}}, {{X, X}}}};
2032
2033 return table[vertex];
2034 }
2036 {
2037 static constexpr ndarray<unsigned int, 6, 2> table = {
2038 {{{0, 0}}, {{0, 1}}, {{0, 2}}, {{0, 3}}, {{1, 2}}, {{X, X}}}};
2039
2040 return table[vertex];
2041 }
2043 {
2044 static constexpr ndarray<unsigned int, 6, 2> table = {
2045 {{{0, 1}}, {{0, 0}}, {{0, 2}}, {{1, 0}}, {{1, 1}}, {{1, 2}}}};
2046
2047 return table[vertex];
2048 }
2050 {
2052 vertex);
2053 }
2054 default:
2056 }
2057
2058 return {};
2059}
2060
2061
2062
2063inline std::array<unsigned int, 2>
2065 const unsigned int line) const
2066{
2067 AssertIndexRange(line, n_lines());
2068
2069 switch (this->kind)
2070 {
2075 {
2077 break;
2078 }
2080 {
2081 static const std::array<unsigned int, 2> table[6] = {
2082 {{0, 0}}, {{0, 1}}, {{0, 2}}, {{1, 1}}, {{1, 2}}, {{2, 1}}};
2083
2084 return table[line];
2085 }
2087 {
2088 static const std::array<unsigned int, 2> table[8] = {{{0, 0}},
2089 {{0, 1}},
2090 {{0, 2}},
2091 {{0, 3}},
2092 {{1, 2}},
2093 {{2, 1}},
2094 {{1, 1}},
2095 {{2, 2}}};
2096
2097 return table[line];
2098 }
2100 {
2101 static const std::array<unsigned int, 2> table[9] = {{{0, 0}},
2102 {{0, 2}},
2103 {{0, 1}},
2104 {{1, 0}},
2105 {{1, 1}},
2106 {{1, 2}},
2107 {{2, 0}},
2108 {{2, 1}},
2109 {{3, 1}}};
2110
2111 return table[line];
2112 }
2114 {
2116 }
2117 default:
2119 }
2120
2121 return {};
2122}
2123
2124
2125inline unsigned int
2127 const unsigned int vertex) const
2128{
2130 AssertIndexRange(line, n_lines());
2131
2132 switch (this->kind)
2133 {
2136 return vertex;
2138 {
2139 static constexpr ndarray<unsigned int, 3, 2> table = {
2140 {{{0, 1}}, {{1, 2}}, {{2, 0}}}};
2141 return table[line][vertex];
2142 }
2144 {
2145 static constexpr ndarray<unsigned int, 4, 2> table = {
2146 {{{0, 2}}, {{1, 3}}, {{0, 1}}, {{2, 3}}}};
2147 return table[line][vertex];
2148 }
2150 {
2151 static constexpr ndarray<unsigned int, 6, 2> table = {
2152 {{{0, 1}}, {{1, 2}}, {{2, 0}}, {{0, 3}}, {{1, 3}}, {{2, 3}}}};
2153 return table[line][vertex];
2154 }
2156 {
2157 static constexpr ndarray<unsigned int, 8, 2> table = {{{{0, 2}},
2158 {{1, 3}},
2159 {{0, 1}},
2160 {{2, 3}},
2161 {{4, 0}},
2162 {{1, 4}},
2163 {{2, 4}},
2164 {{4, 3}}}};
2165 return table[line][vertex];
2166 }
2168 {
2169 static constexpr ndarray<unsigned int, 9, 2> table = {{{{1, 0}},
2170 {{2, 1}},
2171 {{0, 2}},
2172 {{3, 4}},
2173 {{4, 5}},
2174 {{5, 3}},
2175 {{0, 3}},
2176 {{1, 4}},
2177 {{2, 5}}}};
2178 return table[line][vertex];
2179 }
2181 {
2182 // first four lines comprise the bottom face, next four are the top,
2183 // and the last four are 'bottom to top'
2184 static constexpr ndarray<unsigned int, 12, 2> table = {{{{0, 2}},
2185 {{1, 3}},
2186 {{0, 1}},
2187 {{2, 3}},
2188 {{4, 6}},
2189 {{5, 7}},
2190 {{4, 5}},
2191 {{6, 7}},
2192 {{0, 4}},
2193 {{1, 5}},
2194 {{2, 6}},
2195 {{3, 7}}}};
2196 return table[line][vertex];
2197 }
2198
2199 default:
2201 }
2202
2204}
2205
2206
2207inline unsigned int
2209 const unsigned int face,
2210 const unsigned int line,
2211 const unsigned char combined_face_orientation) const
2212{
2213 AssertIndexRange(face, n_faces());
2215
2216 static constexpr unsigned int X = numbers::invalid_unsigned_int;
2217
2218 switch (this->kind)
2219 {
2221 {
2223 break;
2224 }
2226 {
2227 const auto [face_orientation, face_rotation, face_flip] =
2228 internal::split_face_orientation(combined_face_orientation);
2229
2231 face, line, face_orientation, face_flip, face_rotation);
2232 }
2234 {
2235 return face;
2236 }
2238 {
2239 const auto [face_orientation, face_rotation, face_flip] =
2240 internal::split_face_orientation(combined_face_orientation);
2241
2243 face, line, face_orientation, face_flip, face_rotation);
2244 }
2246 {
2247 static constexpr ndarray<unsigned int, 4, 3> table = {
2248 {{{0, 1, 2}}, {{0, 3, 4}}, {{2, 5, 3}}, {{1, 4, 5}}}};
2249
2250 return table[face][standard_to_real_face_line(
2251 line, face, combined_face_orientation)];
2252 }
2254 {
2255 static constexpr ndarray<unsigned int, 5, 4> table = {
2256 {{{0, 1, 2, 3}},
2257 {{0, 6, 4, X}},
2258 {{1, 5, 7, X}},
2259 {{2, 4, 5, X}},
2260 {{3, 7, 6, X}}}};
2261
2262 return table[face][standard_to_real_face_line(
2263 line, face, combined_face_orientation)];
2264 }
2266 {
2267 static constexpr ndarray<unsigned int, 5, 4> table = {
2268 {{{0, 2, 1, X}},
2269 {{3, 4, 5, X}},
2270 {{6, 7, 0, 3}},
2271 {{7, 8, 1, 4}},
2272 {{8, 6, 5, 2}}}};
2273
2274 return table[face][standard_to_real_face_line(
2275 line, face, combined_face_orientation)];
2276 }
2278 {
2279 const auto [face_orientation, face_rotation, face_flip] =
2280 internal::split_face_orientation(combined_face_orientation);
2281
2283 face, line, face_orientation, face_flip, face_rotation);
2284 }
2285 default:
2287 }
2288
2290}
2291
2292
2293
2294inline unsigned int
2296 const unsigned int face,
2297 const unsigned int vertex,
2298 const unsigned char combined_face_orientation) const
2299{
2300 AssertIndexRange(face, n_faces());
2302 // TODO: once the default orientation is switched to 0 then we can remove this
2303 // special case for 1D.
2304 if (get_dimension() == 1)
2305 Assert(combined_face_orientation ==
2307 ExcMessage("In 1D, all faces must have the default orientation."));
2308 else
2309 AssertIndexRange(combined_face_orientation, n_face_orientations(face));
2310
2311 switch (this->kind)
2312 {
2314 {
2316 break;
2317 }
2319 {
2320 const auto [face_orientation, face_rotation, face_flip] =
2321 internal::split_face_orientation(combined_face_orientation);
2322
2324 face, vertex, face_orientation, face_flip, face_rotation);
2325 }
2327 {
2328 static constexpr ndarray<unsigned int, 3, 2> table = {
2329 {{{0, 1}}, {{1, 2}}, {{2, 0}}}};
2330
2331 return table[face][combined_face_orientation ==
2333 vertex :
2334 (1 - vertex)];
2335 }
2337 {
2338 const auto [face_orientation, face_rotation, face_flip] =
2339 internal::split_face_orientation(combined_face_orientation);
2340
2342 face, vertex, face_orientation, face_flip, face_rotation);
2343 }
2345 {
2346 static constexpr ndarray<unsigned int, 4, 3> table = {
2347 {{{0, 1, 2}}, {{1, 0, 3}}, {{0, 2, 3}}, {{2, 1, 3}}}};
2348
2349 return table[face][standard_to_real_face_vertex(
2350 vertex, face, combined_face_orientation)];
2351 }
2353 {
2354 constexpr auto X = numbers::invalid_unsigned_int;
2355 static constexpr ndarray<unsigned int, 5, 4> table = {
2356 {{{0, 1, 2, 3}},
2357 {{0, 2, 4, X}},
2358 {{3, 1, 4, X}},
2359 {{1, 0, 4, X}},
2360 {{2, 3, 4, X}}}};
2361
2362 return table[face][standard_to_real_face_vertex(
2363 vertex, face, combined_face_orientation)];
2364 }
2366 {
2367 constexpr auto X = numbers::invalid_unsigned_int;
2368 static constexpr ndarray<unsigned int, 6, 4> table = {
2369 {{{1, 0, 2, X}},
2370 {{3, 4, 5, X}},
2371 {{0, 1, 3, 4}},
2372 {{1, 2, 4, 5}},
2373 {{2, 0, 5, 3}}}};
2374
2375 return table[face][standard_to_real_face_vertex(
2376 vertex, face, combined_face_orientation)];
2377 }
2379 {
2380 const auto [face_orientation, face_rotation, face_flip] =
2381 internal::split_face_orientation(combined_face_orientation);
2382
2384 face, vertex, face_orientation, face_flip, face_rotation);
2385 }
2386 default:
2388 }
2389
2391}
2392
2393
2394
2395template <int dim>
2398 const unsigned int vertex) const
2399{
2400 return this->template vertex<dim>(
2402}
2403
2404
2405
2406inline unsigned int
2408 const unsigned int vertex,
2409 const unsigned int face,
2410 const unsigned char face_orientation) const
2411{
2412 AssertIndexRange(face, n_faces());
2414
2415 switch (this->kind)
2416 {
2419 break;
2421 Assert(face_orientation == default_combined_face_orientation(),
2422 ExcMessage(
2423 "In 1D, all faces must have the default orientation."));
2424 return vertex;
2427 return line_vertex_permutations[face_orientation][vertex];
2429 return triangle_vertex_permutations[face_orientation][vertex];
2431 // face 0 is a quadrilateral
2432 if (face == 0)
2433 return quadrilateral_vertex_permutations[face_orientation][vertex];
2434 else
2435 return triangle_vertex_permutations[face_orientation][vertex];
2437 // faces 0 and 1 are triangles
2438 if (face > 1)
2439 return quadrilateral_vertex_permutations[face_orientation][vertex];
2440 else
2441 return triangle_vertex_permutations[face_orientation][vertex];
2443 return quadrilateral_vertex_permutations[face_orientation][vertex];
2444 default:
2446 }
2447
2450}
2451
2452
2453
2454inline unsigned int
2456 const unsigned int line,
2457 const unsigned int face,
2458 const unsigned char combined_face_orientation) const
2459{
2460 AssertIndexRange(face, n_faces());
2462
2463 static constexpr ndarray<unsigned int, 6, 3> triangle_table = {{{{2, 1, 0}},
2464 {{0, 1, 2}},
2465 {{1, 0, 2}},
2466 {{2, 0, 1}},
2467 {{0, 2, 1}},
2468 {{1, 2, 0}}}};
2469
2470
2471 switch (this->kind)
2472 {
2478 break;
2480 return triangle_table[combined_face_orientation][line];
2482 if (face == 0) // The quadrilateral face
2483 {
2484 const auto [face_orientation, face_rotation, face_flip] =
2485 internal::split_face_orientation(combined_face_orientation);
2486
2488 face_orientation,
2489 face_flip,
2490 face_rotation);
2491 }
2492 else // One of the triangular faces
2493 {
2494 return triangle_table[combined_face_orientation][line];
2495 }
2497 if (face > 1) // One of the quadrilateral faces
2498 {
2499 const auto [face_orientation, face_rotation, face_flip] =
2500 internal::split_face_orientation(combined_face_orientation);
2501
2503 face_orientation,
2504 face_flip,
2505 face_rotation);
2506 }
2507 else // One of the triangular faces
2508 return triangle_table[combined_face_orientation][line];
2510 {
2511 static constexpr ndarray<unsigned int, 8, 4> table = {
2512 {{{2, 3, 0, 1}},
2513 {{0, 1, 2, 3}},
2514 {{0, 1, 3, 2}},
2515 {{3, 2, 0, 1}},
2516 {{3, 2, 1, 0}},
2517 {{1, 0, 3, 2}},
2518 {{1, 0, 2, 3}},
2519 {{2, 3, 1, 0}}}};
2520 return table[combined_face_orientation][line];
2521 }
2522 default:
2524 }
2525
2527}
2528
2529
2530
2531namespace ReferenceCells
2532{
2533 template <int dim>
2534 inline constexpr const ReferenceCell &
2536 {
2537 switch (dim)
2538 {
2539 case 0:
2541 case 1:
2542 return ReferenceCells::Line;
2543 case 2:
2545 case 3:
2547 default:
2549 }
2551 }
2552
2553
2554
2555 template <int dim>
2556 inline constexpr const ReferenceCell &
2558 {
2559 switch (dim)
2560 {
2561 case 0:
2563 case 1:
2564 return ReferenceCells::Line;
2565 case 2:
2567 case 3:
2569 default:
2571 }
2573 }
2574} // namespace ReferenceCells
2575
2576
2577inline ReferenceCell
2578ReferenceCell::n_vertices_to_type(const int dim, const unsigned int n_vertices)
2579{
2580 AssertIndexRange(dim, 4);
2582
2583 const auto X = ReferenceCells::Invalid;
2584 static const ndarray<ReferenceCell, 4, 9> table = {
2585 {// dim 0
2586 {{X, ReferenceCells::Vertex, X, X, X, X, X, X, X}},
2587 // dim 1
2588 {{X, X, ReferenceCells::Line, X, X, X, X, X, X}},
2589 // dim 2
2590 {{X,
2591 X,
2592 X,
2595 X,
2596 X,
2597 X,
2598 X}},
2599 // dim 3
2600 {{X,
2601 X,
2602 X,
2603 X,
2607 X,
2610 ExcMessage("The combination of dim = " + std::to_string(dim) +
2611 " and n_vertices = " + std::to_string(n_vertices) +
2612 " does not correspond to a known reference cell type."));
2613 return table[dim][n_vertices];
2614}
2615
2616
2617
2618template <int dim>
2619inline double
2621 const unsigned int i) const
2622{
2625 switch (this->kind)
2626 {
2632 // see also BarycentricPolynomials<2>::compute_value
2634 {
2635 switch (i)
2636 {
2637 case 0:
2638 return 1.0 - xi[std::min(0, dim - 1)] -
2639 xi[std::min(1, dim - 1)];
2640 case 1:
2641 return xi[std::min(0, dim - 1)];
2642 case 2:
2643 return xi[std::min(1, dim - 1)];
2644 default:
2646 }
2647 }
2648 // see also BarycentricPolynomials<3>::compute_value
2650 {
2651 switch (i)
2652 {
2653 case 0:
2654 return 1.0 - xi[std::min(0, dim - 1)] -
2655 xi[std::min(1, dim - 1)] - xi[std::min(2, dim - 1)];
2656 case 1:
2657 return xi[std::min(0, dim - 1)];
2658 case 2:
2659 return xi[std::min(1, dim - 1)];
2660 case 3:
2661 return xi[std::min(2, dim - 1)];
2662 default:
2664 }
2665 }
2666 // see also ScalarLagrangePolynomialPyramid::compute_value()
2668 {
2669 const double Q14 = 0.25;
2670
2671 const double r = xi[std::min(0, dim - 1)];
2672 const double s = xi[std::min(1, dim - 1)];
2673 const double t = xi[std::min(2, dim - 1)];
2674
2675 const double ratio =
2676 (std::fabs(t - 1.0) > 1.0e-14 ? (r * s * t) / (1.0 - t) : 0.0);
2677
2678 if (i == 0)
2679 return Q14 * ((1.0 - r) * (1.0 - s) - t + ratio);
2680 if (i == 1)
2681 return Q14 * ((1.0 + r) * (1.0 - s) - t - ratio);
2682 if (i == 2)
2683 return Q14 * ((1.0 - r) * (1.0 + s) - t - ratio);
2684 if (i == 3)
2685 return Q14 * ((1.0 + r) * (1.0 + s) - t + ratio);
2686 else
2687 return t;
2688 }
2689 // see also ScalarLagrangePolynomialWedge::compute_value()
2692 .d_linear_shape_function<2>(Point<2>(xi[std::min(0, dim - 1)],
2693 xi[std::min(1, dim - 1)]),
2694 i % 3) *
2696 .d_linear_shape_function<1>(Point<1>(xi[std::min(2, dim - 1)]),
2697 i / 3);
2698 default:
2700 }
2701
2702 return 0.0;
2703}
2704
2705
2706
2707template <int dim>
2708inline Tensor<1, dim>
2710 const unsigned int i) const
2711{
2713 switch (this->kind)
2714 {
2720 // see also BarycentricPolynomials<2>::compute_grad()
2722 switch (i)
2723 {
2724 case 0:
2725 return Point<dim>(-1.0, -1.0);
2726 case 1:
2727 return Point<dim>(+1.0, +0.0);
2728 case 2:
2729 return Point<dim>(+0.0, +1.0);
2730 default:
2732 }
2733 default:
2735 }
2736
2737 return Point<dim>(+0.0, +0.0, +0.0);
2738}
2739
2740
2741
2742inline double
2744{
2745 switch (this->kind)
2746 {
2748 return 0;
2750 return 1;
2752 return 1. / 2.;
2754 return 1;
2756 return 1. / 6.;
2758 return 4. / 3.;
2760 return 1. / 2.;
2762 return 1;
2763 default:
2765 }
2766
2767 return 0.0;
2768}
2769
2770
2771
2772template <int dim>
2773inline Point<dim>
2775{
2777
2778 switch (this->kind)
2779 {
2781 return Point<dim>();
2783 return Point<dim>(1. / 2.);
2785 return Point<dim>(1. / 3., 1. / 3.);
2787 return Point<dim>(1. / 2., 1. / 2.);
2789 return Point<dim>(1. / 4., 1. / 4., 1. / 4.);
2791 return Point<dim>(0, 0, 1. / 4.);
2793 return Point<dim>(1. / 3, 1. / 3, 1. / 2.);
2795 return Point<dim>(1. / 2., 1. / 2., 1. / 2.);
2796 default:
2798 }
2799
2800 return Point<dim>();
2801}
2802
2803
2804
2805template <int dim>
2806inline bool
2807ReferenceCell::contains_point(const Point<dim> &p, const double tolerance) const
2808{
2810
2811 // Introduce abbreviations to silence compiler warnings about invalid
2812 // array accesses (that can't happen, of course, but the compiler
2813 // doesn't know that).
2814 constexpr unsigned int x_coordinate = 0;
2815 constexpr unsigned int y_coordinate = (dim >= 2 ? 1 : 0);
2816 constexpr unsigned int z_coordinate = (dim >= 3 ? 2 : 0);
2817
2818 switch (this->kind)
2819 {
2821 {
2822 // Vertices are special cases in that they do not actually
2823 // have coordinates. Error out if this function is called
2824 // with a vertex:
2825 Assert(false,
2826 ExcMessage("Vertices are zero-dimensional objects and "
2827 "as a consequence have no coordinates. You "
2828 "cannot meaningfully ask whether a point is "
2829 "inside a vertex (within a certain tolerance) "
2830 "without coordinate values."));
2831 return false;
2832 }
2836 {
2837 for (unsigned int d = 0; d < dim; ++d)
2838 if ((p[d] < -tolerance) || (p[d] > 1 + tolerance))
2839 return false;
2840 return true;
2841 }
2844 {
2845 // First make sure that we are in the first quadrant or octant
2846 for (unsigned int d = 0; d < dim; ++d)
2847 if (p[d] < -tolerance)
2848 return false;
2849
2850 // Now we also need to make sure that we are below the diagonal line
2851 // or plane that delineates the simplex. This diagonal is given by
2852 // sum(p[d])<=1, and a diagonal a distance eps away is given by
2853 // sum(p[d])<=1+eps*sqrt(d). (For example, the point at (1,1) is a
2854 // distance of 1/sqrt(2) away from the diagonal. That is, its
2855 // sum satisfies
2856 // sum(p[d]) = 2 <= 1 + (1/sqrt(2)) * sqrt(2)
2857 // in other words, it satisfies the predicate with eps=1/sqrt(2).)
2858 double sum = 0;
2859 for (unsigned int d = 0; d < dim; ++d)
2860 sum += p[d];
2861 return (sum <= 1 + tolerance * std::sqrt(1. * dim));
2862 }
2864 {
2865 // A pyramid only lives in the upper half-space:
2866 if (p[z_coordinate] < -tolerance)
2867 return false;
2868
2869 // It also only lives in the space below z=1:
2870 if (p[z_coordinate] > 1 + tolerance)
2871 return false;
2872
2873 // Within what's left of the space, a pyramid is a cone that tapers
2874 // towards the top. First compute the distance of the point to the
2875 // axis in the max norm (this is the right norm because the vertices
2876 // of the pyramid are at points +/-1, +/-1):
2877 const double distance_from_axis =
2878 std::max(std::fabs(p[x_coordinate]), std::fabs(p[y_coordinate]));
2879
2880 // We are inside the pyramid if the distance from the axis is less
2881 // than (1-z)
2882 return (distance_from_axis <= 1 + tolerance - p[z_coordinate]);
2883 }
2885 {
2886 // The wedge we use is a triangle extruded into the third
2887 // dimension by one unit. So we can use the same logic as for
2888 // triangles above (i.e., for the simplex above, using dim==2)
2889 // and then check the third dimension separately.
2890
2891 if ((p[x_coordinate] < -tolerance) || (p[y_coordinate] < -tolerance))
2892 return false;
2893
2894 const double sum = p[x_coordinate] + p[y_coordinate];
2895 if (sum > 1 + tolerance * std::sqrt(2.0))
2896 return false;
2897
2898 if (p[z_coordinate] < -tolerance)
2899 return false;
2900 if (p[z_coordinate] > 1 + tolerance)
2901 return false;
2902
2903 return true;
2904 }
2905 default:
2907 }
2908
2909 return false;
2910}
2911
2912
2913
2914template <int dim>
2915inline Tensor<1, dim>
2917 const unsigned int i) const
2918{
2920 AssertIndexRange(i, dim - 1);
2921
2922 switch (this->kind)
2923 {
2931 {
2932 AssertIndexRange(face_no, 3);
2933 static const std::array<Tensor<1, dim>, 3> table = {
2934 {Point<dim>(1, 0),
2935 Point<dim>(-std::sqrt(0.5), +std::sqrt(0.5)),
2936 Point<dim>(0, -1)}};
2937
2938 return table[face_no];
2939 }
2941 {
2942 AssertIndexRange(face_no, 4);
2943 static const ndarray<Tensor<1, dim>, 4, 2> table = {
2944 {{{Point<dim>(0, 1, 0), Point<dim>(1, 0, 0)}},
2945 {{Point<dim>(1, 0, 0), Point<dim>(0, 0, 1)}},
2946 {{Point<dim>(0, 0, 1), Point<dim>(0, 1, 0)}},
2947 {{Point<dim>(-std::pow(1.0 / 3.0, 1.0 / 4.0),
2948 +std::pow(1.0 / 3.0, 1.0 / 4.0),
2949 0),
2950 Point<dim>(-std::pow(1.0 / 3.0, 1.0 / 4.0),
2951 0,
2952 +std::pow(1.0 / 3.0, 1.0 / 4.0))}}}};
2953
2954 return table[face_no][i];
2955 }
2957 {
2958 AssertIndexRange(face_no, 5);
2959 static const ndarray<Tensor<1, dim>, 5, 2> table = {
2960 {{{Point<dim>(0, 1, 0), Point<dim>(1, 0, 0)}},
2961 {{Point<dim>(+1.0 / sqrt(2.0), 0, +1.0 / sqrt(2.0)),
2962 Point<dim>(0, 1, 0)}},
2963 {{Point<dim>(+1.0 / sqrt(2.0), 0, -1.0 / sqrt(2.0)),
2964 Point<dim>(0, 1, 0)}},
2965 {{Point<dim>(1, 0, 0),
2966 Point<dim>(0, +1.0 / sqrt(2.0), +1.0 / sqrt(2.0))}},
2967 {{Point<dim>(1, 0, 0),
2968 Point<dim>(0, +1.0 / sqrt(2.0), -1.0 / sqrt(2.0))}}}};
2969
2970 return table[face_no][i];
2971 }
2973 {
2974 AssertIndexRange(face_no, 5);
2975 static const ndarray<Tensor<1, dim>, 5, 2> table = {
2976 {{{Point<dim>(0, 1, 0), Point<dim>(1, 0, 0)}},
2977 {{Point<dim>(1, 0, 0), Point<dim>(0, 1, 0)}},
2978 {{Point<dim>(1, 0, 0), Point<dim>(0, 0, 1)}},
2979 {{Point<dim>(-1 / std::sqrt(2.0), +1 / std::sqrt(2.0), 0),
2980 Point<dim>(0, 0, 1)}},
2981 {{Point<dim>(0, 0, 1), Point<dim>(0, 1, 0)}}}};
2982
2983 return table[face_no][i];
2984 }
2985 default:
2987 }
2988
2989
2990 return {};
2991}
2992
2993
2994
2995template <int dim>
2996inline Tensor<1, dim>
2997ReferenceCell::unit_normal_vectors(const unsigned int face_no) const
2998{
2999 AssertDimension(dim, this->get_dimension());
3000
3001 if (is_hyper_cube())
3002 {
3005 }
3006 else if (dim == 2)
3007 {
3009
3010 // Return the rotated vector
3011 return cross_product_2d(unit_tangential_vectors<dim>(face_no, 0));
3012 }
3013 else if (dim == 3)
3014 {
3015 return cross_product_3d(unit_tangential_vectors<dim>(face_no, 0),
3016 unit_tangential_vectors<dim>(face_no, 1));
3017 }
3018
3020
3021 return {};
3022}
3023
3024
3025
3026inline unsigned int
3027ReferenceCell::n_face_orientations(const unsigned int face_no) const
3028{
3029 AssertIndexRange(face_no, n_faces());
3030 if (get_dimension() == 1)
3031 return 1;
3032 if (get_dimension() == 2)
3033 return 2;
3035 return 8;
3036 else if (face_reference_cell(face_no) == ReferenceCells::Triangle)
3037 return 6;
3038
3041}
3042
3043
3044
3045inline bool
3047 const unsigned int line,
3048 const unsigned int face,
3049 const unsigned char combined_face_orientation,
3050 const bool line_orientation) const
3051{
3052 if (*this == ReferenceCells::Hexahedron)
3053 {
3054 static constexpr ::ndarray<bool, 2, 8> bool_table{
3055 {{{true, true, false, true, false, false, true, false}},
3056 {{true, true, true, false, false, false, false, true}}}};
3057
3058 return (line_orientation ==
3059 bool_table[line / 2][combined_face_orientation]);
3060 }
3061 else if (*this == ReferenceCells::Tetrahedron)
3062 {
3063 static constexpr unsigned int X = numbers::invalid_unsigned_int;
3064 static constexpr ::ndarray<unsigned int, 4, 3> combined_lines{
3065 {{{0, 0, 0}}, {{X, 0, 1}}, {{X, 0, X}}, {{X, X, X}}}};
3066
3067 const auto combined_line = combined_lines[face][line];
3068
3069 Assert(combined_line != X,
3070 ExcMessage(
3071 "This function can only be called for following face-line "
3072 "combinations: (0,0), (0,1), (0,2), (1,1), (1,2), (2,1),"));
3073
3074 static constexpr ::ndarray<bool, 2, 6> bool_table{
3075 {{{false, true, false, true, false, true}},
3076 {{true, false, true, false, true, false}}}};
3077
3078 return (line_orientation ==
3079 bool_table[combined_line][combined_face_orientation]);
3080 }
3081 else
3082 // TODO: This might actually be wrong for some of the other
3083 // kinds of objects. We should check this
3084 return true;
3085}
3086
3087
3088
3089namespace internal
3090{
3091 template <typename T>
3093 {
3094 public:
3108
3112 virtual ~NoPermutation() noexcept override = default;
3113
3117 virtual void
3118 print_info(std::ostream &out) const override
3119 {
3120 out << '[';
3121
3122 const unsigned int n_vertices = entity_type.n_vertices();
3123
3124 for (unsigned int i = 0; i < n_vertices; ++i)
3125 {
3126 out << vertices_0[i];
3127 if (i + 1 != n_vertices)
3128 out << ',';
3129 }
3130
3131 out << "] is not a valid permutation of [";
3132
3133 for (unsigned int i = 0; i < n_vertices; ++i)
3134 {
3135 out << vertices_1[i];
3136 if (i + 1 != n_vertices)
3137 out << ',';
3138 }
3139
3140 out << "]." << std::endl;
3141 }
3142
3147
3152
3157 };
3158
3170 << "The reference-cell type used on this cell (" << arg1.to_string()
3171 << ") does not match the reference-cell type of the finite element "
3172 << "associated with this cell (" << arg2.to_string() << "). "
3173 << "Did you accidentally use simplex elements on hypercube meshes "
3174 << "(or the other way around), or are you using a mixed mesh and "
3175 << "assigned a simplex element to a hypercube cell (or the other "
3176 << "way around) via the active_fe_index?");
3177} // namespace internal
3178
3179
3180
3181template <typename T, std::size_t N>
3182inline unsigned char
3183ReferenceCell::compute_orientation(const std::array<T, N> &vertices_0,
3184 const std::array<T, N> &vertices_1) const
3185{
3186 Assert(N >= n_vertices(),
3187 ExcMessage("The number of array elements must be equal to or "
3188 "greater than the number of vertices of the cell "
3189 "referenced by this object."));
3190
3191 // Call the non-deprecated function, taking care of calling it only with
3192 // those array elements that we actually care about (see the note
3193 // in the documentation about the arguments potentially being
3194 // larger arrays than necessary).
3196 make_array_view(vertices_0.begin(), vertices_0.begin() + n_vertices()),
3197 make_array_view(vertices_1.begin(), vertices_1.begin() + n_vertices()));
3198}
3199
3200
3201
3202template <typename T>
3203unsigned char
3205 const ArrayView<const T> &vertices_0,
3206 const ArrayView<const T> &vertices_1) const
3207{
3208 Assert(vertices_0.size() == n_vertices(),
3209 ExcMessage("The number of array elements must be equal to "
3210 "the number of vertices of the cell "
3211 "referenced by this object."));
3212 Assert(vertices_1.size() == n_vertices(),
3213 ExcMessage("The number of array elements must be equal to "
3214 "the number of vertices of the cell "
3215 "referenced by this object."));
3216
3217 auto compute_orientation = [&](const auto &table) {
3218 for (unsigned char o = 0; o < table.size(); ++o)
3219 {
3220 bool match = true;
3221 for (unsigned int j = 0; j < table[o].size(); ++j)
3222 match = (match && vertices_0[j] == vertices_1[table[o][j]]);
3223
3224 if (match)
3225 return o;
3226 }
3227
3228 Assert(false, (internal::NoPermutation<T>(*this, vertices_0, vertices_1)));
3229 return std::numeric_limits<unsigned char>::max();
3230 };
3231
3232 switch (this->kind)
3233 {
3235 // TODO: we can get rid of this special-case and use
3236 // vertex_vertex_permutations once we make 0 the default orientation.
3238 if (vertices_0[0] == vertices_1[0])
3240 break;
3247 default:
3249 }
3250
3251 Assert(false, (internal::NoPermutation<T>(*this, vertices_0, vertices_1)));
3252 return std::numeric_limits<unsigned char>::max();
3253}
3254
3255
3256
3257template <typename T, std::size_t N>
3258inline std::array<T, N>
3260 const std::array<T, N> &vertices,
3261 const unsigned int orientation) const
3262{
3263 Assert(N >= n_vertices(),
3264 ExcMessage("The number of array elements must be equal to or "
3265 "greater than the number of vertices of the cell "
3266 "referenced by this object."));
3267
3268 // Call the non-deprecated function, taking care of calling it only with
3269 // those array elements that we actually care about (see the note
3270 // in the documentation about the arguments potentially being
3271 // larger arrays than necessary).
3272 const auto permutation = permute_by_combined_orientation(
3273 make_array_view(vertices.begin(), vertices.begin() + n_vertices()),
3274 orientation);
3275
3276 std::array<T, N> temp;
3277 std::copy(permutation.begin(), permutation.end(), temp.begin());
3278
3279 return temp;
3280}
3281
3282
3283
3284template <typename T>
3285boost::container::small_vector<T, 8>
3288 const unsigned char orientation) const
3289{
3291 boost::container::small_vector<T, 8> result(n_vertices());
3292
3293 auto permute = [&](const auto &table) {
3294 AssertIndexRange(orientation, table.size());
3295 for (unsigned int j = 0; j < table[orientation].size(); ++j)
3296 result[j] = vertices[table[orientation][j]];
3297 };
3298
3299 switch (this->kind)
3300 {
3302 // TODO: we can get rid of this special-case and use
3303 // vertex_vertex_permutations once we make 0 the default orientation.
3304 std::copy(vertices.begin(), vertices.end(), result.begin());
3305 break;
3307 permute(line_vertex_permutations);
3308 break;
3311 break;
3314 break;
3315 default:
3317 }
3318
3319 return result;
3320}
3321
3322
3323
3324inline unsigned char
3326 const unsigned char orientation) const
3327{
3328 switch (this->kind)
3329 {
3331 // Things are always default-oriented in 1D
3332 return orientation;
3333
3335 // the 1d orientations are the identity and a flip: i.e., the identity
3336 // and an involutory mapping
3337 return orientation;
3338
3340 {
3341 AssertIndexRange(orientation, 6);
3342 constexpr std::array<unsigned char, 6> inverses{{0, 1, 2, 5, 4, 3}};
3343 return inverses[orientation];
3344 }
3346 {
3347 AssertIndexRange(orientation, 8);
3348 constexpr std::array<unsigned char, 8> inverses{
3349 {0, 1, 2, 7, 4, 5, 6, 3}};
3350 return inverses[orientation];
3351 }
3352 default:
3354 }
3355
3356 return std::numeric_limits<unsigned char>::max();
3357}
3358
3359
3360
3361template <>
3362unsigned int
3364 const std::array<unsigned, 0> &node_indices,
3365 const std::array<unsigned, 0> &nodes_per_direction,
3366 const bool legacy_format) const;
3367
3368template <>
3369unsigned int
3371 const std::array<unsigned, 1> &node_indices,
3372 const std::array<unsigned, 1> &nodes_per_direction,
3373 const bool legacy_format) const;
3374
3375template <>
3376unsigned int
3378 const std::array<unsigned, 2> &node_indices,
3379 const std::array<unsigned, 2> &nodes_per_direction,
3380 const bool legacy_format) const;
3381
3382template <>
3383unsigned int
3385 const std::array<unsigned, 3> &node_indices,
3386 const std::array<unsigned, 3> &nodes_per_direction,
3387 const bool legacy_format) const;
3388
3390
3391#endif
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:949
std::size_t size() const
Definition array_view.h:684
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
static constexpr Point< dim, Number > unit_vector(const unsigned int i)
ArrayView< const unsigned int > faces_for_given_vertex(const unsigned int vertex_index) const
unsigned int child_cell_on_face(const unsigned int face, const unsigned int subface, const unsigned char face_orientation=default_combined_face_orientation()) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
unsigned int n_vertices() const
boost::container::small_vector< T, 8 > permute_by_combined_orientation(const ArrayView< const T > &vertices, const unsigned char orientation) const
Quadrature< dim > get_gauss_type_quadrature(const unsigned n_points_1d) const
bool is_hyper_cube() const
unsigned int vtk_linear_type() const
void serialize(Archive &archive, const unsigned int)
unsigned char get_inverse_combined_orientation(const unsigned char orientation) const
unsigned int n_isotropic_refinement_choices() const
unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const unsigned char face_orientation) const
static constexpr unsigned char default_combined_face_orientation()
unsigned char compute_orientation(const std::array< T, N > &vertices_0, const std::array< T, N > &vertices_1) const
bool standard_vs_true_line_orientation(const unsigned int line, const unsigned int face, const unsigned char face_orientation, const bool line_orientation) const
std::array< unsigned int, 2 > standard_vertex_to_face_and_vertex_index(const unsigned int vertex) const
Point< dim > vertex(const unsigned int v) const
unsigned int standard_to_real_face_vertex(const unsigned int vertex, const unsigned int face, const unsigned char face_orientation) const
const Quadrature< dim > & get_nodal_type_quadrature() const
Point< dim > face_vertex_location(const unsigned int face, const unsigned int vertex) const
constexpr ::ndarray< unsigned int, 12, 4 > new_isotropic_child_face_lines(const unsigned int refinement_choice) const
static constexpr ndarray< unsigned int, 1, 1 > vertex_vertex_permutations
std_cxx20::ranges::iota_view< unsigned int, unsigned int > isotropic_child_indices() const
std::array< unsigned int, 2 > standard_line_to_face_and_line_index(const unsigned int line) const
constexpr ::ndarray< unsigned int, 12, 4, 2 > new_isotropic_child_face_line_vertices(const unsigned int refinement_choice) 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
static constexpr ndarray< unsigned int, 2, 2 > line_vertex_permutations
unsigned int vtk_lexicographic_to_node_index(const std::array< unsigned, dim > &node_indices, const std::array< unsigned, dim > &nodes_per_direction, const bool legacy_format) const
std::array< T, N > permute_according_orientation(const std::array< T, N > &vertices, const unsigned int orientation) const
constexpr ::ndarray< unsigned int, 8, 4 > new_isotropic_child_cell_vertices(const unsigned int refinement_choice) const
bool contains_point(const Point< dim > &p, const double tolerance=0) const
constexpr ::ndarray< unsigned int, 8, 6 > new_isotropic_child_cell_faces(const unsigned int refinement_choice) const
unsigned int n_isotropic_children() const
unsigned int gmsh_element_type() const
unsigned int n_face_orientations(const unsigned int face_no) const
static constexpr ndarray< unsigned int, 8, 4 > quadrilateral_vertex_permutations
double volume() const
unsigned int n_faces() const
unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex) 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
unsigned int vtk_vertex_to_deal_vertex(const unsigned int vertex_index) const
IsotropicRefinementChoice get_isotropic_refinement_choice(const unsigned int ref_choice) const
std::unique_ptr< Mapping< dim, spacedim > > get_default_mapping(const unsigned int degree) const
Quadrature< dim > get_midpoint_quadrature() const
unsigned int vtk_quadratic_type() const
static constexpr unsigned char reversed_combined_line_orientation()
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
unsigned int n_children(const RefinementCase< dim > ref_case=RefinementCase< dim >::isotropic_refinement) 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
static constexpr std::size_t memory_consumption()
friend std::istream & operator>>(std::istream &in, ReferenceCell &reference_cell)
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
unsigned char get_combined_orientation(const ArrayView< const T > &vertices_0, const ArrayView< const T > &vertices_1) const
constexpr bool operator==(const ReferenceCell &type) const
Point< dim > barycenter() 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
friend std::ostream & operator<<(std::ostream &out, const ReferenceCell &reference_cell)
Point< dim > closest_point(const Point< dim > &p) const
static constexpr ndarray< unsigned int, 6, 3 > triangle_vertex_permutations
Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i) const
virtual void print_info(std::ostream &out) const override
const ReferenceCell entity_type
virtual ~NoPermutation() noexcept override=default
const ArrayView< const T > vertices_0
const ArrayView< const T > vertices_1
NoPermutation(const ReferenceCell &entity_type, const ArrayView< const T > &vertices_0, const ArrayView< const T > &vertices_1)
#define DEAL_II_DEPRECATED
Definition config.h:207
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
Point< 3 > vertices[4]
static ::ExceptionBase & ExcNonMatchingReferenceCellTypes(ReferenceCell arg1, ReferenceCell arg2)
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
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()
constexpr ReferenceCell make_reference_cell_from_int(const std::uint8_t kind)
std::tuple< bool, bool, bool > split_face_orientation(const unsigned char combined_face_orientation)
static const unsigned int invalid_unsigned_int
Definition types.h:220
boost::integer_range< IncrementableType > iota_view
Definition iota_view.h:45
STL namespace.
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(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:107
IsotropicRefinementChoice
std::istream & operator>>(std::istream &in, ReferenceCell &reference_cell)
std::ostream & operator<<(std::ostream &out, const ReferenceCell &reference_cell)
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 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 n_children(const RefinementCase< dim > &refinement_case)
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)