16 #ifndef dealii__geometry_info_h 17 #define dealii__geometry_info_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/point.h> 24 DEAL_II_NAMESPACE_OPEN
78 operator unsigned int ()
const;
131 isotropic_refinement =
static_cast<unsigned char>(-1)
176 isotropic_refinement = cut_x
222 cut_xy = cut_x | cut_y,
224 isotropic_refinement = cut_xy
270 cut_xy = cut_x | cut_y,
272 cut_xz = cut_x | cut_z,
273 cut_yz = cut_y | cut_z,
274 cut_xyz = cut_x | cut_y | cut_z,
276 isotropic_refinement = cut_xyz
325 operator unsigned char ()
const;
360 static std::size_t memory_consumption ();
366 template <
class Archive>
367 void serialize(Archive &ar,
368 const unsigned int version);
375 <<
"The refinement flags given (" << arg1 <<
") contain set bits that do not " 376 <<
"make sense for the space dimension of the object to which they are applied.");
383 unsigned char value :
422 case_isotropic =
static_cast<unsigned char>(-1)
448 case_isotropic = case_none
476 case_isotropic = case_none
507 case_isotropic = case_x
627 case_isotropic = case_xy
663 operator unsigned char ()
const;
668 static std::size_t memory_consumption ();
675 <<
"The subface case given (" << arg1 <<
") does not make sense " 676 <<
"for the space dimension of the object to which they are applied.");
683 unsigned char value :
724 static const unsigned int max_children_per_cell = 1;
729 static const unsigned int faces_per_cell = 0;
738 static const unsigned int max_children_per_face = 0;
750 static const unsigned int vertices_per_cell = 1;
758 static const unsigned int vertices_per_face = 0;
763 static const unsigned int lines_per_face = 0;
768 static const unsigned int quads_per_face = 0;
773 static const unsigned int lines_per_cell = 0;
778 static const unsigned int quads_per_cell = 0;
783 static const unsigned int hexes_per_cell = 0;
1325 static const unsigned int max_children_per_cell = 1 << dim;
1330 static const unsigned int faces_per_cell = 2 * dim;
1339 static const unsigned int max_children_per_face =
GeometryInfo<dim-1>::max_children_per_cell;
1344 static const unsigned int vertices_per_cell = 1 << dim;
1349 static const unsigned int vertices_per_face =
GeometryInfo<dim-1>::vertices_per_cell;
1354 static const unsigned int lines_per_face
1360 static const unsigned int quads_per_face
1372 static const unsigned int lines_per_cell
1383 static const unsigned int quads_per_cell
1390 static const unsigned int hexes_per_cell
1411 static const unsigned int ucd_to_deal[vertices_per_cell];
1426 static const unsigned int dx_to_deal[vertices_per_cell];
1438 static const unsigned int vertex_to_face[vertices_per_cell][dim];
1468 const unsigned int subface_no);
1478 const unsigned int face_no,
1479 const bool face_orientation =
true,
1480 const bool face_flip =
false,
1481 const bool face_rotation =
false);
1490 min_cell_refinement_case_for_face_refinement
1492 const unsigned int face_no,
1493 const bool face_orientation =
true,
1494 const bool face_flip =
false,
1495 const bool face_rotation =
false);
1504 const unsigned int line_no);
1512 min_cell_refinement_case_for_line_refinement(
const unsigned int line_no);
1563 const unsigned int face,
1564 const unsigned int subface,
1565 const bool face_orientation =
true,
1566 const bool face_flip =
false,
1567 const bool face_rotation =
false,
1586 line_to_cell_vertices (
const unsigned int line,
1587 const unsigned int vertex);
1611 face_to_cell_vertices (
const unsigned int face,
1612 const unsigned int vertex,
1613 const bool face_orientation =
true,
1614 const bool face_flip =
false,
1615 const bool face_rotation =
false);
1630 face_to_cell_lines (
const unsigned int face,
1631 const unsigned int line,
1632 const bool face_orientation =
true,
1633 const bool face_flip =
false,
1634 const bool face_rotation =
false);
1647 standard_to_real_face_vertex (
const unsigned int vertex,
1648 const bool face_orientation =
true,
1649 const bool face_flip =
false,
1650 const bool face_rotation =
false);
1663 real_to_standard_face_vertex (
const unsigned int vertex,
1664 const bool face_orientation =
true,
1665 const bool face_flip =
false,
1666 const bool face_rotation =
false);
1679 standard_to_real_face_line (
const unsigned int line,
1680 const bool face_orientation =
true,
1681 const bool face_flip =
false,
1682 const bool face_rotation =
false);
1695 real_to_standard_face_line (
const unsigned int line,
1696 const bool face_orientation =
true,
1697 const bool face_flip =
false,
1698 const bool face_rotation =
false);
1707 unit_cell_vertex (
const unsigned int vertex);
1731 cell_to_child_coordinates (
const Point<dim> &p,
1732 const unsigned int child_index,
1743 child_to_cell_coordinates (
const Point<dim> &p,
1744 const unsigned int child_index,
1795 d_linear_shape_function (
const Point<dim> &xi,
1796 const unsigned int i);
1804 d_linear_shape_function_gradient (
const Point<dim> &xi,
1805 const unsigned int i);
1858 template <
int spacedim>
1860 alternating_form_at_vertices
1861 #ifndef DEAL_II_CONSTEXPR_BUG 1866 Tensor<spacedim-dim,spacedim> *forms);
1878 static const unsigned int unit_normal_direction[faces_per_cell];
1896 static const int unit_normal_orientation[faces_per_cell];
1903 static const unsigned int opposite_face[faces_per_cell];
1911 <<
"The coordinates must satisfy 0 <= x_i <= 1, " 1912 <<
"but here we have x_i=" << arg1);
1919 <<
"RefinementCase<dim> " << arg1 <<
": face " << arg2
1920 <<
" has no subface " << arg3);
1931 #ifndef DEAL_II_MEMBER_ARRAY_SPECIALIZATION_BUG 1965 const unsigned int i);
1970 const unsigned int i);
1975 const unsigned int i);
1994 object (static_cast<Object>(object_dimension))
1999 GeometryPrimitive::operator
unsigned int ()
const 2001 return static_cast<unsigned int>(object);
2011 SubfaceCase<dim>::SubfaceCase (
const typename SubfacePossibilities<dim>::Possibilities subface_possibility)
2013 value (subface_possibility)
2019 SubfaceCase<dim>::operator
unsigned char ()
const 2033 Assert (
false, ExcInternalError());
2034 return static_cast<unsigned char>(-1);
2043 const unsigned int dim = 1;
2044 Assert (i < dim, ExcIndexRange(i, 0, dim));
2057 const unsigned int dim = 2;
2058 Assert (i < dim, ExcIndexRange(i, 0, dim));
2073 const unsigned int dim = 3;
2074 Assert (i < dim, ExcIndexRange(i, 0, dim));
2099 value (refinement_case)
2107 ExcInvalidRefinementCase (refinement_case));
2116 value (refinement_case)
2124 ExcInvalidRefinementCase (refinement_case));
2181 template <
class Archive>
2187 unsigned char uchar_value = value;
2189 value = uchar_value;
2203 return Point<1>(
static_cast<double>(vertex));
2216 return Point<2>(vertex%2, vertex/2);
2229 return Point<3>(vertex%2, vertex/2%2, vertex/4);
2239 Assert(
false, ExcNotImplemented());
2251 Assert ((p[0] >= 0) && (p[0] <= 1), ExcInvalidCoordinate(p[0]));
2253 return (p[0] <= 0.5 ? 0 : 1);
2263 Assert ((p[0] >= 0) && (p[0] <= 1), ExcInvalidCoordinate(p[0]));
2264 Assert ((p[1] >= 0) && (p[1] <= 1), ExcInvalidCoordinate(p[1]));
2266 return (p[0] <= 0.5 ?
2267 (p[1] <= 0.5 ? 0 : 2) :
2268 (p[1] <= 0.5 ? 1 : 3));
2278 Assert ((p[0] >= 0) && (p[0] <= 1), ExcInvalidCoordinate(p[0]));
2279 Assert ((p[1] >= 0) && (p[1] <= 1), ExcInvalidCoordinate(p[1]));
2280 Assert ((p[2] >= 0) && (p[2] <= 1), ExcInvalidCoordinate(p[2]));
2282 return (p[0] <= 0.5 ?
2284 (p[2] <= 0.5 ? 0 : 4) :
2285 (p[2] <= 0.5 ? 2 : 6)) :
2287 (p[2] <= 0.5 ? 1 : 5) :
2288 (p[2] <= 0.5 ? 3 : 7)));
2297 Assert(
false, ExcNotImplemented());
2308 const unsigned int child_index,
2313 ExcIndexRange (child_index, 0, 2));
2315 ExcInternalError());
2327 const unsigned int child_index,
2335 switch (refine_case)
2352 Assert(
false, ExcInternalError());
2364 const unsigned int child_index,
2377 switch (refine_case)
2397 if (child_index%2==1)
2399 if (child_index/2==1)
2409 if (child_index/2==1)
2411 if (child_index%2==1)
2417 if (child_index%2==1)
2419 if (child_index/2==1)
2427 Assert(
false, ExcInternalError());
2439 const unsigned int ,
2453 const unsigned int child_index,
2458 ExcIndexRange (child_index, 0, 2));
2460 ExcInternalError());
2472 const unsigned int child_index,
2485 switch (refine_case)
2503 if (child_index%2==1)
2505 if (child_index/2==1)
2515 if (child_index/2==1)
2517 if (child_index%2==1)
2523 if (child_index%2==1)
2525 if (child_index/2==1)
2535 Assert(
false, ExcInternalError());
2547 const unsigned int child_index,
2554 switch (refine_case)
2571 Assert(
false, ExcInternalError());
2583 const unsigned int ,
2597 return (p[0] >= 0.) && (p[0] <= 1.);
2607 return (p[0] >= 0.) && (p[0] <= 1.) &&
2608 (p[1] >= 0.) && (p[1] <= 1.);
2618 return (p[0] >= 0.) && (p[0] <= 1.) &&
2619 (p[1] >= 0.) && (p[1] <= 1.) &&
2620 (p[2] >= 0.) && (p[2] <= 1.);
2629 return (p[0] >= -eps) && (p[0] <= 1.+eps);
2640 const double l = -eps, u = 1+eps;
2641 return (p[0] >= l) && (p[0] <= u) &&
2642 (p[1] >= l) && (p[1] <= u);
2653 const double l = -eps, u = 1.0+eps;
2654 return (p[0] >= l) && (p[0] <= u) &&
2655 (p[1] >= l) && (p[1] <= u) &&
2656 (p[2] >= l) && (p[2] <= u);
2661 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 child_cell_from_point(const Point< dim > &p)
static bool is_inside_unit_cell(const Point< dim > &p)
GeometryPrimitive(const Object object)
void serialize(Archive &ar, const unsigned int version)
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)
#define AssertThrow(cond, exc)
static const unsigned int vertices_per_cell
RefinementCase operator~() const
RefinementCase operator|(const RefinementCase &r) const
#define DeclException1(Exception1, type1, outsequence)
#define Assert(cond, exc)
UpdateFlags operator|(UpdateFlags f1, UpdateFlags f2)
RefinementCase operator&(const RefinementCase &r) const
static std::size_t memory_consumption()
UpdateFlags operator&(UpdateFlags f1, UpdateFlags f2)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static Tensor< 1, dim > d_linear_shape_function_gradient(const Point< dim > &xi, const unsigned int i)
static RefinementCase cut_axis(const unsigned int i)