Reference documentation for deal.II version 8.4.2
geometry_info.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2015 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__geometry_info_h
17 #define dealii__geometry_info_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/point.h>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 
45 {
46 public:
53  enum Object
54  {
55  vertex = 0,
56  line = 1,
57  quad = 2,
58  hex = 3
59  };
60 
65  GeometryPrimitive (const Object object);
66 
72  GeometryPrimitive (const unsigned int object_dimension);
73 
78  operator unsigned int () const;
79 
80 private:
85 };
86 
87 
88 
89 
103 template <int dim>
105 {
128  {
129  no_refinement= 0,
130 
131  isotropic_refinement = static_cast<unsigned char>(-1)
132  };
133 };
134 
135 
136 
147 template <>
149 {
172  {
173  no_refinement= 0,
174  cut_x = 1,
175 
176  isotropic_refinement = cut_x
177  };
178 };
179 
180 
181 
193 template <>
195 {
218  {
219  no_refinement= 0,
220  cut_x = 1,
221  cut_y = 2,
222  cut_xy = cut_x | cut_y,
223 
224  isotropic_refinement = cut_xy
225  };
226 };
227 
228 
229 
241 template <>
243 {
266  {
267  no_refinement= 0,
268  cut_x = 1,
269  cut_y = 2,
270  cut_xy = cut_x | cut_y,
271  cut_z = 4,
272  cut_xz = cut_x | cut_z,
273  cut_yz = cut_y | cut_z,
274  cut_xyz = cut_x | cut_y | cut_z,
275 
276  isotropic_refinement = cut_xyz
277  };
278 };
279 
280 
281 
292 template <int dim>
294 {
295 public:
299  RefinementCase ();
300 
305  RefinementCase (const typename RefinementPossibilities<dim>::Possibilities refinement_case);
306 
312  explicit RefinementCase (const unsigned char refinement_case);
313 
325  operator unsigned char () const;
326 
331  RefinementCase operator | (const RefinementCase &r) const;
332 
337  RefinementCase operator & (const RefinementCase &r) const;
338 
346  RefinementCase operator ~ () const;
347 
348 
354  static
355  RefinementCase cut_axis (const unsigned int i);
356 
360  static std::size_t memory_consumption ();
361 
366  template <class Archive>
367  void serialize(Archive &ar,
368  const unsigned int version);
369 
373  DeclException1 (ExcInvalidRefinementCase,
374  int,
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.");
377 
378 private:
383 unsigned char value :
384  (dim > 0 ? dim : 1);
385 };
386 
387 
388 namespace internal
389 {
390 
391 
412  template <int dim>
414  {
419  {
420  case_none = 0,
421 
422  case_isotropic = static_cast<unsigned char>(-1)
423  };
424  };
425 
426 
436  template <>
438  {
445  {
446  case_none = 0,
447 
448  case_isotropic = case_none
449  };
450  };
451 
452 
453 
464  template <>
466  {
473  {
474  case_none = 0,
475 
476  case_isotropic = case_none
477  };
478  };
479 
480 
481 
493  template <>
495  {
503  {
504  case_none = 0,
505  case_x = 1,
506 
507  case_isotropic = case_x
508  };
509  };
510 
511 
512 
605  template <>
607  {
615  {
616  case_none = 0,
617  case_x = 1,
618  case_x1y = 2,
619  case_x2y = 3,
620  case_x1y2y = 4,
621  case_y = 5,
622  case_y1x = 6,
623  case_y2x = 7,
624  case_y1x2x = 8,
625  case_xy = 9,
626 
627  case_isotropic = case_xy
628  };
629  };
630 
631 
632 
633 
641  template <int dim>
642  class SubfaceCase : public SubfacePossibilities<dim>
643  {
644  public:
650  SubfaceCase (const typename SubfacePossibilities<dim>::Possibilities subface_possibility);
651 
663  operator unsigned char () const;
664 
668  static std::size_t memory_consumption ();
669 
673  DeclException1 (ExcInvalidSubfaceCase,
674  int,
675  << "The subface case given (" << arg1 << ") does not make sense "
676  << "for the space dimension of the object to which they are applied.");
677 
678  private:
683  unsigned char value :
684  (dim == 3 ? 4 : 1);
685  };
686 
687 } // namespace internal
688 
689 
690 
691 template <int dim> struct GeometryInfo;
692 
693 
694 
695 
713 template <>
714 struct GeometryInfo<0>
715 {
716 
724  static const unsigned int max_children_per_cell = 1;
725 
729  static const unsigned int faces_per_cell = 0;
730 
738  static const unsigned int max_children_per_face = 0;
739 
745  static unsigned int n_children(const RefinementCase<0> &refinement_case);
746 
750  static const unsigned int vertices_per_cell = 1;
751 
758  static const unsigned int vertices_per_face = 0;
759 
763  static const unsigned int lines_per_face = 0;
764 
768  static const unsigned int quads_per_face = 0;
769 
773  static const unsigned int lines_per_cell = 0;
774 
778  static const unsigned int quads_per_cell = 0;
779 
783  static const unsigned int hexes_per_cell = 0;
784 };
785 
786 
787 
788 
789 
1314 template <int dim>
1315 struct GeometryInfo
1316 {
1317 
1325  static const unsigned int max_children_per_cell = 1 << dim;
1326 
1330  static const unsigned int faces_per_cell = 2 * dim;
1331 
1339  static const unsigned int max_children_per_face = GeometryInfo<dim-1>::max_children_per_cell;
1340 
1344  static const unsigned int vertices_per_cell = 1 << dim;
1345 
1349  static const unsigned int vertices_per_face = GeometryInfo<dim-1>::vertices_per_cell;
1350 
1354  static const unsigned int lines_per_face
1355  = GeometryInfo<dim-1>::lines_per_cell;
1356 
1360  static const unsigned int quads_per_face
1361  = GeometryInfo<dim-1>::quads_per_cell;
1362 
1372  static const unsigned int lines_per_cell
1373  = (2*GeometryInfo<dim-1>::lines_per_cell +
1374  GeometryInfo<dim-1>::vertices_per_cell);
1375 
1383  static const unsigned int quads_per_cell
1384  = (2*GeometryInfo<dim-1>::quads_per_cell +
1385  GeometryInfo<dim-1>::lines_per_cell);
1386 
1390  static const unsigned int hexes_per_cell
1391  = (2*GeometryInfo<dim-1>::hexes_per_cell +
1392  GeometryInfo<dim-1>::quads_per_cell);
1393 
1411  static const unsigned int ucd_to_deal[vertices_per_cell];
1412 
1426  static const unsigned int dx_to_deal[vertices_per_cell];
1427 
1438  static const unsigned int vertex_to_face[vertices_per_cell][dim];
1439 
1444  static
1445  unsigned int
1446  n_children(const RefinementCase<dim> &refinement_case);
1447 
1452  static
1453  unsigned int
1454  n_subfaces(const internal::SubfaceCase<dim> &subface_case);
1455 
1465  static
1466  double
1467  subface_ratio(const internal::SubfaceCase<dim> &subface_case,
1468  const unsigned int subface_no);
1469 
1475  static
1476  RefinementCase<dim-1>
1477  face_refinement_case (const RefinementCase<dim> &cell_refinement_case,
1478  const unsigned int face_no,
1479  const bool face_orientation = true,
1480  const bool face_flip = false,
1481  const bool face_rotation = false);
1482 
1488  static
1490  min_cell_refinement_case_for_face_refinement
1491  (const RefinementCase<dim-1> &face_refinement_case,
1492  const unsigned int face_no,
1493  const bool face_orientation = true,
1494  const bool face_flip = false,
1495  const bool face_rotation = false);
1496 
1501  static
1503  line_refinement_case(const RefinementCase<dim> &cell_refinement_case,
1504  const unsigned int line_no);
1505 
1510  static
1512  min_cell_refinement_case_for_line_refinement(const unsigned int line_no);
1513 
1560  static
1561  unsigned int
1562  child_cell_on_face (const RefinementCase<dim> &ref_case,
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,
1568  const RefinementCase<dim-1> &face_refinement_case
1570 
1584  static
1585  unsigned int
1586  line_to_cell_vertices (const unsigned int line,
1587  const unsigned int vertex);
1588 
1609  static
1610  unsigned int
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);
1616 
1628  static
1629  unsigned int
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);
1635 
1645  static
1646  unsigned int
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);
1651 
1661  static
1662  unsigned int
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);
1667 
1677  static
1678  unsigned int
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);
1683 
1693  static
1694  unsigned int
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);
1699 
1705  static
1706  Point<dim>
1707  unit_cell_vertex (const unsigned int vertex);
1708 
1718  static
1719  unsigned int
1720  child_cell_from_point (const Point<dim> &p);
1721 
1729  static
1730  Point<dim>
1731  cell_to_child_coordinates (const Point<dim> &p,
1732  const unsigned int child_index,
1733  const RefinementCase<dim> refine_case
1735 
1741  static
1742  Point<dim>
1743  child_to_cell_coordinates (const Point<dim> &p,
1744  const unsigned int child_index,
1745  const RefinementCase<dim> refine_case
1747 
1752  static
1753  bool
1754  is_inside_unit_cell (const Point<dim> &p);
1755 
1767  static
1768  bool
1769  is_inside_unit_cell (const Point<dim> &p,
1770  const double eps);
1771 
1776  static
1777  Point<dim>
1778  project_to_unit_cell (const Point<dim> &p);
1779 
1785  static
1786  double
1787  distance_to_unit_cell (const Point<dim> &p);
1788 
1793  static
1794  double
1795  d_linear_shape_function (const Point<dim> &xi,
1796  const unsigned int i);
1797 
1802  static
1804  d_linear_shape_function_gradient (const Point<dim> &xi,
1805  const unsigned int i);
1806 
1858  template <int spacedim>
1859  static void
1860  alternating_form_at_vertices
1861 #ifndef DEAL_II_CONSTEXPR_BUG
1862  (const Point<spacedim> (&vertices)[vertices_per_cell],
1863  Tensor<spacedim-dim,spacedim> (&forms)[vertices_per_cell]);
1864 #else
1865  (const Point<spacedim> *vertices,
1866  Tensor<spacedim-dim,spacedim> *forms);
1867 #endif
1868 
1878  static const unsigned int unit_normal_direction[faces_per_cell];
1879 
1896  static const int unit_normal_orientation[faces_per_cell];
1897 
1903  static const unsigned int opposite_face[faces_per_cell];
1904 
1905 
1909  DeclException1 (ExcInvalidCoordinate,
1910  double,
1911  << "The coordinates must satisfy 0 <= x_i <= 1, "
1912  << "but here we have x_i=" << arg1);
1913 
1917  DeclException3 (ExcInvalidSubface,
1918  int, int, int,
1919  << "RefinementCase<dim> " << arg1 << ": face " << arg2
1920  << " has no subface " << arg3);
1921 };
1922 
1923 
1924 
1925 
1926 #ifndef DOXYGEN
1927 
1928 
1929 /* -------------- declaration of explicit specializations ------------- */
1930 
1931 #ifndef DEAL_II_MEMBER_ARRAY_SPECIALIZATION_BUG
1932 template <>
1933 const unsigned int GeometryInfo<1>::unit_normal_direction[faces_per_cell];
1934 template <>
1935 const unsigned int GeometryInfo<2>::unit_normal_direction[faces_per_cell];
1936 template <>
1937 const unsigned int GeometryInfo<3>::unit_normal_direction[faces_per_cell];
1938 template <>
1939 const unsigned int GeometryInfo<4>::unit_normal_direction[faces_per_cell];
1940 
1941 template <>
1942 const int GeometryInfo<1>::unit_normal_orientation[faces_per_cell];
1943 template <>
1944 const int GeometryInfo<2>::unit_normal_orientation[faces_per_cell];
1945 template <>
1946 const int GeometryInfo<3>::unit_normal_orientation[faces_per_cell];
1947 template <>
1948 const int GeometryInfo<4>::unit_normal_orientation[faces_per_cell];
1949 
1950 template <>
1951 const unsigned int GeometryInfo<1>::opposite_face[faces_per_cell];
1952 template <>
1953 const unsigned int GeometryInfo<2>::opposite_face[faces_per_cell];
1954 template <>
1955 const unsigned int GeometryInfo<3>::opposite_face[faces_per_cell];
1956 template <>
1957 const unsigned int GeometryInfo<4>::opposite_face[faces_per_cell];
1958 #endif
1959 
1960 
1961 template <>
1965  const unsigned int i);
1966 template <>
1970  const unsigned int i);
1971 template <>
1975  const unsigned int i);
1976 
1977 
1978 
1979 
1980 /* -------------- inline functions ------------- */
1981 
1982 
1983 inline
1985  :
1986  object (object)
1987 {}
1988 
1989 
1990 
1991 inline
1992 GeometryPrimitive::GeometryPrimitive (const unsigned int object_dimension)
1993  :
1994  object (static_cast<Object>(object_dimension))
1995 {}
1996 
1997 
1998 inline
1999 GeometryPrimitive::operator unsigned int () const
2000 {
2001  return static_cast<unsigned int>(object);
2002 }
2003 
2004 
2005 
2006 namespace internal
2007 {
2008 
2009  template <int dim>
2010  inline
2011  SubfaceCase<dim>::SubfaceCase (const typename SubfacePossibilities<dim>::Possibilities subface_possibility)
2012  :
2013  value (subface_possibility)
2014  {}
2015 
2016 
2017  template <int dim>
2018  inline
2019  SubfaceCase<dim>::operator unsigned char () const
2020  {
2021  return value;
2022  }
2023 
2024 
2025 } // namespace internal
2026 
2027 
2028 template <int dim>
2029 inline
2031 RefinementCase<dim>::cut_axis (const unsigned int)
2032 {
2033  Assert (false, ExcInternalError());
2034  return static_cast<unsigned char>(-1);
2035 }
2036 
2037 
2038 template <>
2039 inline
2041 RefinementCase<1>::cut_axis (const unsigned int i)
2042 {
2043  const unsigned int dim = 1;
2044  Assert (i < dim, ExcIndexRange(i, 0, dim));
2045 
2046  static const RefinementCase options[dim] = { RefinementPossibilities<1>::cut_x };
2047  return options[i];
2048 }
2049 
2050 
2051 
2052 template <>
2053 inline
2055 RefinementCase<2>::cut_axis (const unsigned int i)
2056 {
2057  const unsigned int dim = 2;
2058  Assert (i < dim, ExcIndexRange(i, 0, dim));
2059 
2060  static const RefinementCase options[dim] = { RefinementPossibilities<2>::cut_x,
2062  };
2063  return options[i];
2064 }
2065 
2066 
2067 
2068 template <>
2069 inline
2071 RefinementCase<3>::cut_axis (const unsigned int i)
2072 {
2073  const unsigned int dim = 3;
2074  Assert (i < dim, ExcIndexRange(i, 0, dim));
2075 
2076  static const RefinementCase options[dim] = { RefinementPossibilities<3>::cut_x,
2079  };
2080  return options[i];
2081 }
2082 
2083 
2084 
2085 template <int dim>
2086 inline
2088  :
2089  value (RefinementPossibilities<dim>::no_refinement)
2090 {}
2091 
2092 
2093 
2094 template <int dim>
2095 inline
2097 RefinementCase (const typename RefinementPossibilities<dim>::Possibilities refinement_case)
2098  :
2099  value (refinement_case)
2100 {
2101  // check that only those bits of
2102  // the given argument are set that
2103  // make sense for a given space
2104  // dimension
2106  refinement_case,
2107  ExcInvalidRefinementCase (refinement_case));
2108 }
2109 
2110 
2111 
2112 template <int dim>
2113 inline
2114 RefinementCase<dim>::RefinementCase (const unsigned char refinement_case)
2115  :
2116  value (refinement_case)
2117 {
2118  // check that only those bits of
2119  // the given argument are set that
2120  // make sense for a given space
2121  // dimension
2123  refinement_case,
2124  ExcInvalidRefinementCase (refinement_case));
2125 }
2126 
2127 
2128 
2129 template <int dim>
2130 inline
2131 RefinementCase<dim>::operator unsigned char () const
2132 {
2133  return value;
2134 }
2135 
2136 
2137 
2138 template <int dim>
2139 inline
2142 {
2143  return RefinementCase<dim>(static_cast<unsigned char> (value | r.value));
2144 }
2145 
2146 
2147 
2148 template <int dim>
2149 inline
2152 {
2153  return RefinementCase<dim>(static_cast<unsigned char> (value & r.value));
2154 }
2155 
2156 
2157 
2158 template <int dim>
2159 inline
2162 {
2163  return RefinementCase<dim>(static_cast<unsigned char> (
2165 }
2166 
2167 
2168 
2169 
2170 template <int dim>
2171 inline
2172 std::size_t
2174 {
2175  return sizeof(RefinementCase<dim>);
2176 }
2177 
2178 
2179 
2180 template <int dim>
2181 template <class Archive>
2182 void RefinementCase<dim>::serialize (Archive &ar,
2183  const unsigned int)
2184 {
2185  // serialization can't deal with bitfields, so copy from/to a full sized
2186  // unsigned char
2187  unsigned char uchar_value = value;
2188  ar &uchar_value;
2189  value = uchar_value;
2190 }
2191 
2192 
2193 
2194 
2195 template <>
2196 inline
2197 Point<1>
2198 GeometryInfo<1>::unit_cell_vertex (const unsigned int vertex)
2199 {
2200  Assert (vertex < vertices_per_cell,
2201  ExcIndexRange (vertex, 0, vertices_per_cell));
2202 
2203  return Point<1>(static_cast<double>(vertex));
2204 }
2205 
2206 
2207 
2208 template <>
2209 inline
2210 Point<2>
2211 GeometryInfo<2>::unit_cell_vertex (const unsigned int vertex)
2212 {
2213  Assert (vertex < vertices_per_cell,
2214  ExcIndexRange (vertex, 0, vertices_per_cell));
2215 
2216  return Point<2>(vertex%2, vertex/2);
2217 }
2218 
2219 
2220 
2221 template <>
2222 inline
2223 Point<3>
2224 GeometryInfo<3>::unit_cell_vertex (const unsigned int vertex)
2225 {
2226  Assert (vertex < vertices_per_cell,
2227  ExcIndexRange (vertex, 0, vertices_per_cell));
2228 
2229  return Point<3>(vertex%2, vertex/2%2, vertex/4);
2230 }
2231 
2232 
2233 
2234 template <int dim>
2235 inline
2236 Point<dim>
2237 GeometryInfo<dim>::unit_cell_vertex (const unsigned int)
2238 {
2239  Assert(false, ExcNotImplemented());
2240 
2241  return Point<dim> ();
2242 }
2243 
2244 
2245 
2246 template <>
2247 inline
2248 unsigned int
2250 {
2251  Assert ((p[0] >= 0) && (p[0] <= 1), ExcInvalidCoordinate(p[0]));
2252 
2253  return (p[0] <= 0.5 ? 0 : 1);
2254 }
2255 
2256 
2257 
2258 template <>
2259 inline
2260 unsigned int
2262 {
2263  Assert ((p[0] >= 0) && (p[0] <= 1), ExcInvalidCoordinate(p[0]));
2264  Assert ((p[1] >= 0) && (p[1] <= 1), ExcInvalidCoordinate(p[1]));
2265 
2266  return (p[0] <= 0.5 ?
2267  (p[1] <= 0.5 ? 0 : 2) :
2268  (p[1] <= 0.5 ? 1 : 3));
2269 }
2270 
2271 
2272 
2273 template <>
2274 inline
2275 unsigned int
2277 {
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]));
2281 
2282  return (p[0] <= 0.5 ?
2283  (p[1] <= 0.5 ?
2284  (p[2] <= 0.5 ? 0 : 4) :
2285  (p[2] <= 0.5 ? 2 : 6)) :
2286  (p[1] <= 0.5 ?
2287  (p[2] <= 0.5 ? 1 : 5) :
2288  (p[2] <= 0.5 ? 3 : 7)));
2289 }
2290 
2291 
2292 template <int dim>
2293 inline
2294 unsigned int
2296 {
2297  Assert(false, ExcNotImplemented());
2298 
2299  return 0;
2300 }
2301 
2302 
2303 
2304 template <>
2305 inline
2306 Point<1>
2308  const unsigned int child_index,
2309  const RefinementCase<1> refine_case)
2310 
2311 {
2312  Assert (child_index < 2,
2313  ExcIndexRange (child_index, 0, 2));
2314  Assert (refine_case==RefinementCase<1>::cut_x,
2315  ExcInternalError());
2316  (void)refine_case; // removes -Wunused-parameter warning in optimized mode
2317 
2318  return Point<1>(p*2.0-unit_cell_vertex(child_index));
2319 }
2320 
2321 
2322 
2323 template <>
2324 inline
2325 Point<2>
2327  const unsigned int child_index,
2328  const RefinementCase<2> refine_case)
2329 
2330 {
2331  Assert (child_index < GeometryInfo<2>::n_children(refine_case),
2332  ExcIndexRange (child_index, 0, GeometryInfo<2>::n_children(refine_case)));
2333 
2334  Point<2> point=p;
2335  switch (refine_case)
2336  {
2338  point[0]*=2.0;
2339  if (child_index==1)
2340  point[0]-=1.0;
2341  break;
2343  point[1]*=2.0;
2344  if (child_index==1)
2345  point[1]-=1.0;
2346  break;
2348  point*=2.0;
2349  point-=unit_cell_vertex(child_index);
2350  break;
2351  default:
2352  Assert(false, ExcInternalError());
2353  }
2354 
2355  return point;
2356 }
2357 
2358 
2359 
2360 template <>
2361 inline
2362 Point<3>
2364  const unsigned int child_index,
2365  const RefinementCase<3> refine_case)
2366 
2367 {
2368  Assert (child_index < GeometryInfo<3>::n_children(refine_case),
2369  ExcIndexRange (child_index, 0, GeometryInfo<3>::n_children(refine_case)));
2370 
2371  Point<3> point=p;
2372  // there might be a cleverer way to do
2373  // this, but since this function is called
2374  // in very few places for initialization
2375  // purposes only, I don't care at the
2376  // moment
2377  switch (refine_case)
2378  {
2380  point[0]*=2.0;
2381  if (child_index==1)
2382  point[0]-=1.0;
2383  break;
2385  point[1]*=2.0;
2386  if (child_index==1)
2387  point[1]-=1.0;
2388  break;
2390  point[2]*=2.0;
2391  if (child_index==1)
2392  point[2]-=1.0;
2393  break;
2395  point[0]*=2.0;
2396  point[1]*=2.0;
2397  if (child_index%2==1)
2398  point[0]-=1.0;
2399  if (child_index/2==1)
2400  point[1]-=1.0;
2401  break;
2403  // careful, this is slightly
2404  // different from xy and yz due to
2405  // different internal numbering of
2406  // children!
2407  point[0]*=2.0;
2408  point[2]*=2.0;
2409  if (child_index/2==1)
2410  point[0]-=1.0;
2411  if (child_index%2==1)
2412  point[2]-=1.0;
2413  break;
2415  point[1]*=2.0;
2416  point[2]*=2.0;
2417  if (child_index%2==1)
2418  point[1]-=1.0;
2419  if (child_index/2==1)
2420  point[2]-=1.0;
2421  break;
2423  point*=2.0;
2424  point-=unit_cell_vertex(child_index);
2425  break;
2426  default:
2427  Assert(false, ExcInternalError());
2428  }
2429 
2430  return point;
2431 }
2432 
2433 
2434 
2435 template <int dim>
2436 inline
2437 Point<dim>
2439  const unsigned int /*child_index*/,
2440  const RefinementCase<dim> /*refine_case*/)
2441 
2442 {
2443  AssertThrow (false, ExcNotImplemented());
2444  return Point<dim>();
2445 }
2446 
2447 
2448 
2449 template <>
2450 inline
2451 Point<1>
2453  const unsigned int child_index,
2454  const RefinementCase<1> refine_case)
2455 
2456 {
2457  Assert (child_index < 2,
2458  ExcIndexRange (child_index, 0, 2));
2459  Assert (refine_case==RefinementCase<1>::cut_x,
2460  ExcInternalError());
2461  (void)refine_case; // removes -Wunused-parameter warning in optimized mode
2462 
2463  return (p+unit_cell_vertex(child_index))*0.5;
2464 }
2465 
2466 
2467 
2468 template <>
2469 inline
2470 Point<3>
2472  const unsigned int child_index,
2473  const RefinementCase<3> refine_case)
2474 
2475 {
2476  Assert (child_index < GeometryInfo<3>::n_children(refine_case),
2477  ExcIndexRange (child_index, 0, GeometryInfo<3>::n_children(refine_case)));
2478 
2479  Point<3> point=p;
2480  // there might be a cleverer way to do
2481  // this, but since this function is called
2482  // in very few places for initialization
2483  // purposes only, I don't care at the
2484  // moment
2485  switch (refine_case)
2486  {
2488  if (child_index==1)
2489  point[0]+=1.0;
2490  point[0]*=0.5;
2491  break;
2493  if (child_index==1)
2494  point[1]+=1.0;
2495  point[1]*=0.5;
2496  break;
2498  if (child_index==1)
2499  point[2]+=1.0;
2500  point[2]*=0.5;
2501  break;
2503  if (child_index%2==1)
2504  point[0]+=1.0;
2505  if (child_index/2==1)
2506  point[1]+=1.0;
2507  point[0]*=0.5;
2508  point[1]*=0.5;
2509  break;
2511  // careful, this is slightly
2512  // different from xy and yz due to
2513  // different internal numbering of
2514  // children!
2515  if (child_index/2==1)
2516  point[0]+=1.0;
2517  if (child_index%2==1)
2518  point[2]+=1.0;
2519  point[0]*=0.5;
2520  point[2]*=0.5;
2521  break;
2523  if (child_index%2==1)
2524  point[1]+=1.0;
2525  if (child_index/2==1)
2526  point[2]+=1.0;
2527  point[1]*=0.5;
2528  point[2]*=0.5;
2529  break;
2531  point+=unit_cell_vertex(child_index);
2532  point*=0.5;
2533  break;
2534  default:
2535  Assert(false, ExcInternalError());
2536  }
2537 
2538  return point;
2539 }
2540 
2541 
2542 
2543 template <>
2544 inline
2545 Point<2>
2547  const unsigned int child_index,
2548  const RefinementCase<2> refine_case)
2549 {
2550  Assert (child_index < GeometryInfo<2>::n_children(refine_case),
2551  ExcIndexRange (child_index, 0, GeometryInfo<2>::n_children(refine_case)));
2552 
2553  Point<2> point=p;
2554  switch (refine_case)
2555  {
2557  if (child_index==1)
2558  point[0]+=1.0;
2559  point[0]*=0.5;
2560  break;
2562  if (child_index==1)
2563  point[1]+=1.0;
2564  point[1]*=0.5;
2565  break;
2567  point+=unit_cell_vertex(child_index);
2568  point*=0.5;
2569  break;
2570  default:
2571  Assert(false, ExcInternalError());
2572  }
2573 
2574  return point;
2575 }
2576 
2577 
2578 
2579 template <int dim>
2580 inline
2581 Point<dim>
2583  const unsigned int /*child_index*/,
2584  const RefinementCase<dim> /*refine_case*/)
2585 {
2586  AssertThrow (false, ExcNotImplemented());
2587  return Point<dim>();
2588 }
2589 
2590 
2591 
2592 template <>
2593 inline
2594 bool
2596 {
2597  return (p[0] >= 0.) && (p[0] <= 1.);
2598 }
2599 
2600 
2601 
2602 template <>
2603 inline
2604 bool
2606 {
2607  return (p[0] >= 0.) && (p[0] <= 1.) &&
2608  (p[1] >= 0.) && (p[1] <= 1.);
2609 }
2610 
2611 
2612 
2613 template <>
2614 inline
2615 bool
2617 {
2618  return (p[0] >= 0.) && (p[0] <= 1.) &&
2619  (p[1] >= 0.) && (p[1] <= 1.) &&
2620  (p[2] >= 0.) && (p[2] <= 1.);
2621 }
2622 
2623 template <>
2624 inline
2625 bool
2627  const double eps)
2628 {
2629  return (p[0] >= -eps) && (p[0] <= 1.+eps);
2630 }
2631 
2632 
2633 
2634 template <>
2635 inline
2636 bool
2638  const double eps)
2639 {
2640  const double l = -eps, u = 1+eps;
2641  return (p[0] >= l) && (p[0] <= u) &&
2642  (p[1] >= l) && (p[1] <= u);
2643 }
2644 
2645 
2646 
2647 template <>
2648 inline
2649 bool
2651  const double eps)
2652 {
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);
2657 }
2658 
2659 
2660 #endif // DOXYGEN
2661 DEAL_II_NAMESPACE_CLOSE
2662 
2663 #endif
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)
Definition: exceptions.h:358
static const unsigned int vertices_per_cell
RefinementCase operator~() const
RefinementCase operator|(const RefinementCase &r) const
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:542
#define Assert(cond, exc)
Definition: exceptions.h:294
unsigned char value
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)
Definition: exceptions.h:562
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)