16 #ifndef dealii_geometry_info_h 17 #define dealii_geometry_info_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/point.h> 29 DEAL_II_NAMESPACE_OPEN
34 namespace GeometryInfoHelper
42 struct Initializers<1>
44 static constexpr std::array<unsigned int, 2> ucd_to_deal{{0, 1}};
46 static constexpr std::array<unsigned int, 2> unit_normal_direction{
49 static constexpr std::array<int, 2> unit_normal_orientation{{-1, 1}};
51 static constexpr std::array<unsigned int, 2> opposite_face{{1, 0}};
53 static constexpr std::array<unsigned int, 2> dx_to_deal{{0, 1}};
55 static constexpr std::array<std::array<unsigned int, 1>, 2>
56 vertex_to_face{{{{0}}, {{1}}}};
60 struct Initializers<2>
62 static constexpr std::array<unsigned int, 4> ucd_to_deal{{0, 1, 3, 2}};
64 static constexpr std::array<unsigned int, 4> unit_normal_direction{
67 static constexpr std::array<int, 4> unit_normal_orientation{
70 static constexpr std::array<unsigned int, 4> opposite_face{{1, 0, 3, 2}};
72 static constexpr std::array<unsigned int, 4> dx_to_deal{{0, 2, 1, 3}};
74 static constexpr std::array<std::array<unsigned int, 2>, 4>
75 vertex_to_face{{{{0, 2}}, {{1, 2}}, {{0, 3}}, {{1, 3}}}};
79 struct Initializers<3>
81 static constexpr std::array<unsigned int, 8> ucd_to_deal{
82 {0, 1, 5, 4, 2, 3, 7, 6}};
84 static constexpr std::array<unsigned int, 6> unit_normal_direction{
87 static constexpr std::array<int, 6> unit_normal_orientation{
88 {-1, 1, -1, 1, -1, 1}};
90 static constexpr std::array<unsigned int, 6> opposite_face{
93 static constexpr std::array<unsigned int, 8> dx_to_deal{
94 {0, 4, 2, 6, 1, 5, 3, 7}};
96 static constexpr std::array<std::array<unsigned int, 3>, 8>
97 vertex_to_face{{{{0, 2, 4}},
108 struct Initializers<4>
110 static constexpr std::array<unsigned int, 16> ucd_to_deal{
128 static constexpr std::array<unsigned int, 8> unit_normal_direction{
129 {0, 0, 1, 1, 2, 2, 3, 3}};
131 static constexpr std::array<int, 8> unit_normal_orientation{
132 {-1, 1, -1, 1, -1, 1, -1, 1}};
134 static constexpr std::array<unsigned int, 8> opposite_face{
135 {1, 0, 3, 2, 5, 4, 7, 6}};
137 static constexpr std::array<unsigned int, 16> dx_to_deal{
155 static constexpr std::array<std::array<unsigned int, 4>, 16>
289 operator unsigned int()
const;
505 cut_xy = cut_x | cut_y,
581 cut_xy = cut_x | cut_y,
589 cut_xz = cut_x | cut_z,
593 cut_yz = cut_y | cut_z,
597 cut_xyz = cut_x | cut_y | cut_z,
653 operator std::uint8_t()
const;
697 template <
class Archive>
699 serialize(Archive &ar,
const unsigned int version);
707 <<
"The refinement flags given (" << arg1
708 <<
") contain set bits that do not " 709 <<
"make sense for the space dimension of the object to which they are applied.");
716 std::uint8_t
value : (dim > 0 ? dim : 1);
1006 subface_possibility);
1019 operator std::uint8_t()
const;
1024 static constexpr std::size_t
1033 <<
"The subface case given (" << arg1 <<
") does not make sense " 1034 <<
"for the space dimension of the object to which they are applied.");
1041 std::uint8_t
value : (dim == 3 ? 4 : 1);
1174 static const std::array<unsigned int, vertices_per_cell>
dx_to_deal;
1801 static constexpr std::array<unsigned int, vertices_per_cell>
ucd_to_deal =
1802 internal::GeometryInfoHelper::Initializers<dim>::ucd_to_deal;
1817 static constexpr std::array<unsigned int, vertices_per_cell>
dx_to_deal =
1818 internal::GeometryInfoHelper::Initializers<dim>::dx_to_deal;
1832 internal::GeometryInfoHelper::Initializers<dim>::vertex_to_face;
1859 const unsigned int subface_no);
1868 const unsigned int face_no,
1869 const bool face_orientation =
true,
1870 const bool face_flip =
false,
1871 const bool face_rotation =
false);
1881 const unsigned int face_no,
1882 const bool face_orientation =
true,
1883 const bool face_flip =
false,
1884 const bool face_rotation =
false);
1892 const unsigned int line_no);
1949 const unsigned int face,
1950 const unsigned int subface,
1951 const bool face_orientation =
true,
1952 const bool face_flip =
false,
1953 const bool face_rotation =
false,
1995 const unsigned int vertex,
1996 const bool face_orientation =
true,
1997 const bool face_flip =
false,
1998 const bool face_rotation =
false);
2013 const unsigned int line,
2014 const bool face_orientation =
true,
2015 const bool face_flip =
false,
2016 const bool face_rotation =
false);
2029 const bool face_orientation =
true,
2030 const bool face_flip =
false,
2031 const bool face_rotation =
false);
2044 const bool face_orientation =
true,
2045 const bool face_flip =
false,
2046 const bool face_rotation =
false);
2059 const bool face_orientation =
true,
2060 const bool face_flip =
false,
2061 const bool face_rotation =
false);
2074 const bool face_orientation =
true,
2075 const bool face_flip =
false,
2076 const bool face_rotation =
false);
2107 const unsigned int child_index,
2118 const unsigned int child_index,
2223 template <
int spacedim>
2226 #ifndef DEAL_II_CONSTEXPR_BUG 2242 static constexpr std::array<unsigned int, faces_per_cell>
2244 internal::GeometryInfoHelper::Initializers<dim>::unit_normal_direction;
2263 internal::GeometryInfoHelper::Initializers<dim>::unit_normal_orientation;
2271 internal::GeometryInfoHelper::Initializers<dim>::opposite_face;
2279 <<
"The coordinates must satisfy 0 <= x_i <= 1, " 2280 <<
"but here we have x_i=" << arg1);
2289 <<
"RefinementCase<dim> " << arg1 <<
": face " << arg2
2290 <<
" has no subface " << arg3);
2304 const unsigned int i);
2308 const unsigned int i);
2312 const unsigned int i);
2326 : object(static_cast<Object>(object_dimension))
2330 inline GeometryPrimitive::operator
unsigned int()
const 2332 return static_cast<unsigned int>(object);
2342 : value(subface_possibility)
2347 inline SubfaceCase<dim>::operator std::uint8_t()
const 2361 return static_cast<std::uint8_t
>(-1);
2414 : value(refinement_case)
2420 Assert((refinement_case &
2423 ExcInvalidRefinementCase(refinement_case));
2430 : value(refinement_case)
2436 Assert((refinement_case &
2439 ExcInvalidRefinementCase(refinement_case));
2490 template <
class Archive>
2496 std::uint8_t uchar_value = value;
2498 value = uchar_value;
2507 Assert(vertex < vertices_per_cell,
2510 return Point<1>(
static_cast<double>(vertex));
2519 Assert(vertex < vertices_per_cell,
2522 return {
static_cast<double>(vertex % 2), static_cast<double>(vertex / 2)};
2531 Assert(vertex < vertices_per_cell,
2534 return {
static_cast<double>(vertex % 2),
2535 static_cast<double>(vertex / 2 % 2),
2536 static_cast<double>(vertex / 4)};
2556 Assert((p[0] >= 0) && (p[0] <= 1), ExcInvalidCoordinate(p[0]));
2558 return (p[0] <= 0.5 ? 0 : 1);
2567 Assert((p[0] >= 0) && (p[0] <= 1), ExcInvalidCoordinate(p[0]));
2568 Assert((p[1] >= 0) && (p[1] <= 1), ExcInvalidCoordinate(p[1]));
2570 return (p[0] <= 0.5 ? (p[1] <= 0.5 ? 0 : 2) : (p[1] <= 0.5 ? 1 : 3));
2579 Assert((p[0] >= 0) && (p[0] <= 1), ExcInvalidCoordinate(p[0]));
2580 Assert((p[1] >= 0) && (p[1] <= 1), ExcInvalidCoordinate(p[1]));
2581 Assert((p[2] >= 0) && (p[2] <= 1), ExcInvalidCoordinate(p[2]));
2583 return (p[0] <= 0.5 ?
2584 (p[1] <= 0.5 ? (p[2] <= 0.5 ? 0 : 4) : (p[2] <= 0.5 ? 2 : 6)) :
2585 (p[1] <= 0.5 ? (p[2] <= 0.5 ? 1 : 5) : (p[2] <= 0.5 ? 3 : 7)));
2603 const unsigned int child_index,
2611 return Point<1>(p * 2.0 - unit_cell_vertex(child_index));
2619 const unsigned int child_index,
2629 switch (refine_case)
2633 if (child_index == 1)
2638 if (child_index == 1)
2643 point -= unit_cell_vertex(child_index);
2657 const unsigned int child_index,
2672 switch (refine_case)
2676 if (child_index == 1)
2681 if (child_index == 1)
2686 if (child_index == 1)
2692 if (child_index % 2 == 1)
2694 if (child_index / 2 == 1)
2704 if (child_index / 2 == 1)
2706 if (child_index % 2 == 1)
2712 if (child_index % 2 == 1)
2714 if (child_index / 2 == 1)
2719 point -= unit_cell_vertex(child_index);
2734 const unsigned int ,
2747 const unsigned int child_index,
2755 return (p + unit_cell_vertex(child_index)) * 0.5;
2763 const unsigned int child_index,
2778 switch (refine_case)
2781 if (child_index == 1)
2786 if (child_index == 1)
2791 if (child_index == 1)
2796 if (child_index % 2 == 1)
2798 if (child_index / 2 == 1)
2808 if (child_index / 2 == 1)
2810 if (child_index % 2 == 1)
2816 if (child_index % 2 == 1)
2818 if (child_index / 2 == 1)
2824 point += unit_cell_vertex(child_index);
2839 const unsigned int child_index,
2848 switch (refine_case)
2851 if (child_index == 1)
2856 if (child_index == 1)
2861 point += unit_cell_vertex(child_index);
2877 const unsigned int ,
2898 return (p[0] >= 0.) && (p[0] <= 1.);
2907 return (p[0] >= 0.) && (p[0] <= 1.) && (p[1] >= 0.) && (p[1] <= 1.);
2916 return (p[0] >= 0.) && (p[0] <= 1.) && (p[1] >= 0.) && (p[1] <= 1.) &&
2917 (p[2] >= 0.) && (p[2] <= 1.);
2934 return (p[0] >= -eps) && (p[0] <= 1. +
eps);
2943 const double l = -
eps, u = 1 +
eps;
2944 return (p[0] >= l) && (p[0] <= u) && (p[1] >= l) && (p[1] <= u);
2953 const double l = -
eps, u = 1.0 +
eps;
2954 return (p[0] >= l) && (p[0] <= u) && (p[1] >= l) && (p[1] <= u) &&
2955 (p[2] >= l) && (p[2] <= u);
2963 const unsigned int vertex)
2976 const unsigned int vertex)
2978 constexpr
unsigned int cell_vertices[4][2] = {{0, 2}, {1, 3}, {0, 1}, {2, 3}};
2979 return cell_vertices[line][vertex];
2987 const unsigned int vertex)
2992 constexpr
unsigned vertices[lines_per_cell][2] = {{0, 2},
3005 return vertices[line][vertex];
3020 const bool face_orientation,
3021 const bool face_flip,
3022 const bool face_rotation)
3045 constexpr
unsigned int vertex_translation[4][2][2][2] = {
3067 return vertex_translation[vertex][face_orientation][face_flip][face_rotation];
3090 0, 2, 2, 4, 2, 4, 4, 8};
3092 return n_children[ref_case];
3127 0, 2, 3, 3, 4, 2, 3, 3, 4, 4};
3128 return nsubs[subface_case];
3157 switch (subface_case)
3194 const unsigned int subface_no)
3197 switch (subface_case)
3232 if (subface_no == 0)
3279 const unsigned int face_no,
3284 const unsigned int dim = 2;
3305 return ref_cases[cell_refinement_case][face_no / 2];
3313 const unsigned int face_no,
3314 const bool face_orientation,
3316 const bool face_rotation)
3318 const unsigned int dim = 3;
3363 ref_cases[cell_refinement_case][face_no / 2];
3383 return (face_orientation == face_rotation) ? flip[ref_case] : ref_case;
3401 const unsigned int line_no)
3404 const unsigned int dim = 1;
3413 return cell_refinement_case;
3421 const unsigned int line_no)
3424 return face_refinement_case(cell_refinement_case, line_no);
3432 const unsigned int line_no)
3434 const unsigned int dim = 3;
3454 const unsigned int direction[lines_per_cell] = {
3455 1, 1, 0, 0, 1, 1, 0, 0, 2, 2, 2, 2};
3457 return ((cell_refinement_case & cut_one[direction[line_no]]) ?
3487 const unsigned int dim = 1;
3498 const unsigned int face_no,
3503 const unsigned int dim = 2;
3504 Assert(face_refinement_case <
3524 const unsigned int face_no,
3525 const bool face_orientation,
3527 const bool face_rotation)
3529 const unsigned int dim = 3;
3530 Assert(face_refinement_case <
3557 (face_orientation == face_rotation) ? flip[face_refinement_case] :
3558 face_refinement_case;
3577 return face_to_cell[face_no / 2][std_face_ref];
3595 const unsigned int line_no)
3607 const unsigned int line_no)
3609 const unsigned int dim = 2;
3621 const unsigned int line_no)
3623 const unsigned int dim = 3;
3635 return ref_cases[line_no / 2];
3643 const bool face_orientation,
3644 const bool face_flip,
3645 const bool face_rotation)
3668 const unsigned int vertex_translation[4][2][2][2] = {
3690 return vertex_translation[vertex][face_orientation][face_flip][face_rotation];
3713 const bool face_orientation,
3714 const bool face_flip,
3715 const bool face_rotation)
3738 const unsigned int line_translation[4][2][2][2] = {
3760 return line_translation[line][face_orientation][face_flip][face_rotation];
3781 const bool face_orientation,
3782 const bool face_flip,
3783 const bool face_rotation)
3806 const unsigned int line_translation[4][2][2][2] = {
3828 return line_translation[line][face_orientation][face_flip][face_rotation];
3849 const unsigned int face,
3850 const unsigned int subface,
3858 Assert(subface < max_children_per_face,
3869 const unsigned int face,
3870 const unsigned int subface,
3872 const bool face_flip,
3877 Assert(subface < max_children_per_face,
3889 [max_children_per_face] = {
3892 {{0, 0}, {1, 1}, {0, 1}, {0, 1}},
3893 {{0, 1}, {0, 1}, {0, 0}, {1, 1}},
3894 {{0, 2}, {1, 3}, {0, 1}, {2, 3}}
3898 {{0, 0}, {1, 1}, {1, 0}, {1, 0}},
3899 {{1, 0}, {1, 0}, {0, 0}, {1, 1}},
3900 {{2, 0}, {3, 1}, {1, 0}, {3, 2}}
3903 return subcells[face_flip][ref_case - 1][face][subface];
3911 const unsigned int face,
3912 const unsigned int subface,
3913 const bool face_orientation,
3914 const bool face_flip,
3915 const bool face_rotation,
3918 const unsigned int dim = 3;
3959 (face_orientation == face_rotation) ? flip[face_ref_case] : face_ref_case;
3971 const unsigned int subface_exchange[4][2][2][2][4] = {
3974 {{{{0,
e,
e,
e}, {0,
e,
e,
e}}, {{0,
e,
e,
e}, {0,
e,
e,
e}}},
3975 {{{0,
e,
e,
e}, {0,
e,
e,
e}}, {{0,
e,
e,
e}, {0,
e,
e,
e}}}},
3989 {{{{0, 1,
e,
e}, {0, 1,
e,
e}}, {{1, 0,
e,
e}, {1, 0,
e,
e}}},
3990 {{{0, 1,
e,
e}, {0, 1,
e,
e}}, {{1, 0,
e,
e}, {1, 0,
e,
e}}}},
3993 {{{{0, 1,
e,
e}, {1, 0,
e,
e}}, {{1, 0,
e,
e}, {0, 1,
e,
e}}},
3994 {{{0, 1,
e,
e}, {1, 0,
e,
e}}, {{1, 0,
e,
e}, {0, 1,
e,
e}}}},
4016 const unsigned int std_subface =
4017 subface_exchange[face_ref_case][face_orientation][face_flip][face_rotation]
4029 [max_children_per_face] = {
4091 if ((std_face_ref & face_refinement_case(ref_case, face)) ==
4092 face_refinement_case(ref_case, face))
4101 const unsigned int equivalent_iso_subface[4][4] = {
4107 const unsigned int equ_std_subface =
4108 equivalent_iso_subface[std_face_ref][std_subface];
4111 return iso_children[ref_case - 1][face][equ_std_subface];
4118 ExcMessage(
"The face RefineCase is too coarse " 4119 "for the given cell RefineCase."));
4146 const unsigned int line,
4166 const unsigned int line,
4184 const unsigned int line,
4185 const bool face_orientation,
4186 const bool face_flip,
4187 const bool face_rotation)
4192 const unsigned lines[faces_per_cell][lines_per_face] = {
4199 return lines[face][real_to_standard_face_line(
4200 line, face_orientation, face_flip, face_rotation)];
4222 const unsigned int vertex,
4223 const bool face_orientation,
4224 const bool face_flip,
4225 const bool face_rotation)
4242 for (
unsigned int i = 0; i < dim; i++)
4257 double result = 0.0;
4259 for (
unsigned int i = 0; i < dim; i++)
4260 if ((-p[i]) > result)
4262 else if ((p[i] - 1.) > result)
4263 result = (p[i] - 1.);
4273 const unsigned int i)
4282 const double x = xi[0];
4295 const double x = xi[0];
4296 const double y = xi[1];
4300 return (1 - x) * (1 - y);
4313 const double x = xi[0];
4314 const double y = xi[1];
4315 const double z = xi[2];
4319 return (1 - x) * (1 - y) * (1 - z);
4321 return x * (1 - y) * (1 - z);
4323 return (1 - x) * y * (1 - z);
4325 return x * y * (1 - z);
4327 return (1 - x) * (1 - y) * z;
4329 return x * (1 - y) * z;
4331 return (1 - x) * y * z;
4349 const unsigned int i)
4370 const unsigned int i)
4375 const double x = xi[0];
4376 const double y = xi[1];
4380 return Point<2>(-(1 - y), -(1 - x));
4396 const unsigned int i)
4401 const double x = xi[0];
4402 const double y = xi[1];
4403 const double z = xi[2];
4407 return Point<3>(-(1 - y) * (1 - z),
4409 -(1 - x) * (1 - y));
4411 return Point<3>((1 - y) * (1 - z), -x * (1 - z), -x * (1 - y));
4413 return Point<3>(-y * (1 - z), (1 - x) * (1 - z), -(1 - x) * y);
4415 return Point<3>(y * (1 - z), x * (1 - z), -x * y);
4417 return Point<3>(-(1 - y) * z, -(1 - x) * z, (1 - x) * (1 - y));
4419 return Point<3>((1 - y) * z, -x * z, x * (1 - y));
4421 return Point<3>(-y * z, (1 - x) * z, (1 - x) * y);
4423 return Point<3>(y * z, x * z, x * y);
4449 namespace GeometryInfoHelper
4459 result[0] = derivative[0][1];
4460 result[1] = -derivative[0][0];
4471 return cross_product_3d(derivative[0], derivative[1]);
4483 for (
unsigned int i = 0; i < dim; ++i)
4484 jacobian[i] = derivative[i];
4493 template <
int spacedim>
4496 # ifndef DEAL_II_CONSTEXPR_BUG 4532 for (
unsigned int i = 0; i < vertices_per_cell; ++i)
4536 for (
unsigned int j = 0; j < vertices_per_cell; ++j)
4539 d_linear_shape_function_gradient(unit_cell_vertex(i), j);
4540 for (
unsigned int l = 0;
l < dim; ++
l)
4541 derivatives[l] += vertices[j] * grad_phi_j[l];
4544 forms[i] = internal::GeometryInfoHelper::wedge_product(derivatives);
4549 DEAL_II_NAMESPACE_CLOSE
static Point< dim > child_to_cell_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
static unsigned int real_to_standard_face_vertex(const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static const unsigned int invalid_unsigned_int
static constexpr unsigned int vertices_per_face
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)
Number determinant(const SymmetricTensor< 2, dim, Number > &)
static unsigned int child_cell_from_point(const Point< dim > &p)
static constexpr std::array< unsigned int, vertices_per_cell > ucd_to_deal
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 RefinementCase< 1 > line_refinement_case(const RefinementCase< dim > &cell_refinement_case, const unsigned int line_no)
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 constexpr unsigned int max_children_per_face
static constexpr unsigned int max_children_per_cell
static bool is_inside_unit_cell(const Point< dim > &p)
GeometryPrimitive(const Object object)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
void serialize(Archive &ar, const unsigned int version)
static unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex)
static constexpr unsigned int quads_per_face
RefinementCase operator &(const RefinementCase &r) const
static unsigned int real_to_standard_face_line(const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static Point< dim > cell_to_child_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
static Point< dim > unit_cell_vertex(const unsigned int vertex)
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
static constexpr std::size_t memory_consumption()
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcInvalidSubfaceCase(int arg1)
static const std::array< unsigned int, vertices_per_cell > ucd_to_deal
static constexpr unsigned int vertices_per_cell
RefinementCase operator~() const
static constexpr std::array< unsigned int, vertices_per_cell > dx_to_deal
static double distance_to_unit_cell(const Point< dim > &p)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcInvalidSubface(int arg1, int arg2, int arg3)
RefinementCase operator|(const RefinementCase &r) const
static constexpr std::array< unsigned int, faces_per_cell > unit_normal_direction
#define DeclException1(Exception1, type1, outsequence)
static constexpr unsigned int lines_per_cell
#define Assert(cond, exc)
static RefinementCase< dim - 1 > face_refinement_case(const RefinementCase< dim > &cell_refinement_case, const unsigned int face_no, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static constexpr unsigned int faces_per_cell
static constexpr unsigned int hexes_per_cell
static ::ExceptionBase & ExcInvalidCoordinate(double arg1)
static void alternating_form_at_vertices(const Point< spacedim >(&vertices)[vertices_per_cell], Tensor< spacedim - dim, spacedim >(&forms)[vertices_per_cell])
static std::size_t memory_consumption()
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)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
static RefinementCase< dim > min_cell_refinement_case_for_line_refinement(const unsigned int line_no)
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
static Point< dim > project_to_unit_cell(const Point< dim > &p)
SubfaceCase(const typename SubfacePossibilities< dim >::Possibilities subface_possibility)
static constexpr unsigned int quads_per_cell
static const std::array< unsigned int, vertices_per_cell > dx_to_deal
static constexpr std::array< int, faces_per_cell > unit_normal_orientation
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static constexpr unsigned int lines_per_face
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 n_subfaces(const internal::SubfaceCase< dim > &subface_case)
static RefinementCase< dim > min_cell_refinement_case_for_face_refinement(const RefinementCase< dim - 1 > &face_refinement_case, const unsigned int face_no, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static ::ExceptionBase & ExcNotImplemented()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static constexpr std::array< unsigned int, faces_per_cell > opposite_face
static ::ExceptionBase & ExcInvalidRefinementCase(int arg1)
static double subface_ratio(const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
static RefinementCase cut_axis(const unsigned int i)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static constexpr std::array< std::array< unsigned int, dim >, vertices_per_cell > vertex_to_face
static ::ExceptionBase & ExcInternalError()