Reference documentation for deal.II version 9.1.0-pre
tria_accessor.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_tria_accessor_h
17 #define dealii_tria_accessor_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/bounding_box.h>
23 #include <deal.II/base/exceptions.h>
24 #include <deal.II/base/geometry_info.h>
25 #include <deal.II/base/point.h>
26 
27 #include <deal.II/grid/cell_id.h>
28 #include <deal.II/grid/tria_iterator_base.h>
29 #include <deal.II/grid/tria_iterator_selector.h>
30 
31 #include <utility>
32 
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 template <int dim, int spacedim>
37 class Triangulation;
38 template <typename Accessor>
39 class TriaRawIterator;
40 template <typename Accessor>
41 class TriaIterator;
42 template <typename Accessor>
43 class TriaActiveIterator;
44 
45 template <int dim, int spacedim>
46 class Manifold;
47 
48 
49 namespace internal
50 {
51  namespace TriangulationImplementation
52  {
53  template <int dim>
54  class TriaObject;
55  template <typename G>
56  class TriaObjects;
57  struct Implementation;
58  } // namespace TriangulationImplementation
59 
60  namespace TriaAccessorImplementation
61  {
62  struct Implementation;
63 
69  template <int structdim, int dim>
71  {
72  struct type
73  {
77  type() = default;
78 
82  type(const int level)
83  {
84  Assert(level == 0, ExcInternalError());
85  (void)level; // removes -Wunused-parameter warning in optimized mode
86  }
87 
91  operator int() const
92  {
93  return 0;
94  }
95 
96  void
97  operator++() const
98  {
99  Assert(false, ExcInternalError());
100  }
101 
102  void
103  operator--() const
104  {
105  Assert(false, ExcInternalError());
106  }
107  };
108  };
109 
110 
116  template <int dim>
117  struct PresentLevelType<dim, dim>
118  {
119  using type = int;
120  };
121 
122  } // namespace TriaAccessorImplementation
123 } // namespace internal
124 template <int structdim, int dim, int spacedim>
125 class TriaAccessor;
126 template <int dim, int spacedim>
127 class TriaAccessor<0, dim, spacedim>;
128 template <int spacedim>
129 class TriaAccessor<0, 1, spacedim>;
130 
131 // note: the file tria_accessor.templates.h is included at the end of
132 // this file. this includes a lot of templates. originally, this was
133 // only done in debug mode, but led to cyclic reduction problems and
134 // so is now on by default.
135 
136 
141 {
146  "The operation you are attempting can only be performed for "
147  "(cell, face, or edge) iterators that point to valid "
148  "objects. These objects need not necessarily be active, "
149  "i.e., have no children, but they need to be part of a "
150  "triangulation. (The objects pointed to by an iterator "
151  "may -- after coarsening -- also be objects that used "
152  "to be part of a triangulation, but are now no longer "
153  "used. Their memory location may have been retained "
154  "for re-use upon the next mesh refinement, but is "
155  "currently unused.)");
166  "The operation you are attempting can only be performed for "
167  "(cell, face, or edge) iterators that point to 'active' "
168  "objects. 'Active' objects are those that do not have "
169  "children (in the case of cells), or that are part of "
170  "an active cell (in the case of faces or edges). However, "
171  "the object on which you are trying the current "
172  "operation is not 'active' in this sense.");
179  "The operation you are attempting can only be performed for "
180  "(cell, face, or edge) iterators that have children, "
181  "but the object on which you are trying the current "
182  "operation does not have any.");
190  "The operation you are attempting can only be performed for "
191  "(cell, face, or edge) iterators that have a parent object, "
192  "but the object on which you are trying the current "
193  "operation does not have one -- i.e., it is on the "
194  "coarsest level of the triangulation.");
199  int,
200  << "You can only set the child index if the cell does not "
201  << "currently have children registered; or you can clear it. "
202  << "The given index was " << arg1
203  << " (-1 means: clear children).");
207  template <typename AccessorType>
209  AccessorType,
210  << "You tried to dereference an iterator for which this "
211  << "is not possible. More information on this iterator: "
212  << "index=" << arg1.index() << ", state="
213  << (arg1.state() == IteratorState::valid ?
214  "valid" :
215  (arg1.state() == IteratorState::past_the_end ?
216  "past_the_end" :
217  "invalid")));
222  "Iterators can only be compared if they point to the same "
223  "triangulation, or if neither of them are associated "
224  "with a triangulation.");
225  // TODO: Write documentation!
230  // TODO: Write documentation!
250  // TODO: Write documentation!
256  int,
257  << "You can only set the child index of an even numbered child."
258  << "The number of the child given was " << arg1 << ".");
259 } // namespace TriaAccessorExceptions
260 
261 
288 template <int structdim, int dim, int spacedim = dim>
289 class TriaAccessorBase
290 {
291 public:
297  static const unsigned int space_dimension = spacedim;
298 
304  static const unsigned int dimension = dim;
305 
311  static const unsigned int structure_dimension = structdim;
312 
313 protected:
319  using AccessorData = void;
320 
324  TriaAccessorBase(const Triangulation<dim, spacedim> *parent = nullptr,
325  const int level = -1,
326  const int index = -1,
327  const AccessorData * = nullptr);
328 
333 
343  void
344  operator=(const TriaAccessorBase *) = delete;
345 
353  void
354  copy_from(const TriaAccessorBase &);
355 
360  operator=(const TriaAccessorBase &);
361 
368  bool
369  operator<(const TriaAccessorBase &other) const;
370 
371 protected:
375  bool
376  operator==(const TriaAccessorBase &) const;
377 
381  bool
382  operator!=(const TriaAccessorBase &) const;
383 
397  void
398  operator++();
399 
407  void
408  operator--();
418  objects() const;
419 
420 public:
426  using LocalData = void *;
427 
451  int
452  level() const;
453 
480  int
481  index() const;
482 
488  state() const;
489 
495  get_triangulation() const;
496 
500 protected:
505  typename ::internal::TriaAccessorImplementation::
506  PresentLevelType<structdim, dim>::type present_level;
507 
513 
518 
519 private:
520  template <typename Accessor>
521  friend class TriaRawIterator;
522  template <typename Accessor>
523  friend class TriaIterator;
524  template <typename Accessor>
525  friend class TriaActiveIterator;
526 };
527 
528 
529 
551 template <int structdim, int dim, int spacedim = dim>
552 class InvalidAccessor : public TriaAccessorBase<structdim, dim, spacedim>
553 {
554 public:
558  using AccessorData =
560 
568  InvalidAccessor(const Triangulation<dim, spacedim> *parent = nullptr,
569  const int level = -1,
570  const int index = -1,
571  const AccessorData * local_data = nullptr);
572 
581 
586  template <typename OtherAccessor>
587  InvalidAccessor(const OtherAccessor &);
588 
592  void
593  copy_from(const InvalidAccessor &);
594 
598  bool
599  operator==(const InvalidAccessor &) const;
600  bool
601  operator!=(const InvalidAccessor &) const;
602 
606  void
607  operator++() const;
608  void
609  operator--() const;
610 
615  bool
616  used() const;
617 
622  bool
623  has_children() const;
624 
629  manifold_id() const;
630 
635  vertex(const unsigned int i) const;
636 
641  typename ::internal::TriangulationImplementation::
642  Iterators<dim, spacedim>::line_iterator
643  line(const unsigned int i) const;
644 
649  typename ::internal::TriangulationImplementation::
650  Iterators<dim, spacedim>::quad_iterator
651  quad(const unsigned int i) const;
652 };
653 
654 
655 
674 template <int structdim, int dim, int spacedim>
675 class TriaAccessor : public TriaAccessorBase<structdim, dim, spacedim>
676 {
677 public:
681  using AccessorData =
683 
687  TriaAccessor(const Triangulation<dim, spacedim> *parent = nullptr,
688  const int level = -1,
689  const int index = -1,
690  const AccessorData * local_data = nullptr);
691 
704  template <int structdim2, int dim2, int spacedim2>
706 
711  template <int structdim2, int dim2, int spacedim2>
713 
723  void
724  operator=(const TriaAccessor &) = delete;
725 
732  bool
733  used() const;
734 
747  vertex_iterator(const unsigned int i) const;
748 
764  unsigned int
765  vertex_index(const unsigned int i) const;
766 
805  vertex(const unsigned int i) const;
806 
810  typename ::internal::TriangulationImplementation::
811  Iterators<dim, spacedim>::line_iterator
812  line(const unsigned int i) const;
813 
820  unsigned int
821  line_index(const unsigned int i) const;
822 
826  typename ::internal::TriangulationImplementation::
827  Iterators<dim, spacedim>::quad_iterator
828  quad(const unsigned int i) const;
829 
836  unsigned int
837  quad_index(const unsigned int i) const;
860  bool
861  face_orientation(const unsigned int face) const;
862 
872  bool
873  face_flip(const unsigned int face) const;
874 
884  bool
885  face_rotation(const unsigned int face) const;
886 
897  bool
898  line_orientation(const unsigned int line) const;
913  bool
914  has_children() const;
915 
920  unsigned int
921  n_children() const;
922 
936  unsigned int
937  number_of_children() const;
938 
952  unsigned int
953  max_refinement_depth() const;
954 
959  child(const unsigned int i) const;
960 
970  isotropic_child(const unsigned int i) const;
971 
976  refinement_case() const;
977 
983  int
984  child_index(const unsigned int i) const;
985 
991  int
992  isotropic_child_index(const unsigned int i) const;
1015  boundary_id() const;
1016 
1046  void
1047  set_boundary_id(const types::boundary_id) const;
1048 
1079  void
1080  set_all_boundary_ids(const types::boundary_id) const;
1081 
1089  bool
1090  at_boundary() const;
1091 
1101  const Manifold<dim, spacedim> &
1102  get_manifold() const;
1103 
1125  manifold_id() const;
1126 
1144  void
1145  set_manifold_id(const types::manifold_id) const;
1146 
1160  void
1161  set_all_manifold_ids(const types::manifold_id) const;
1162 
1179  bool
1180  user_flag_set() const;
1181 
1187  void
1188  set_user_flag() const;
1189 
1195  void
1196  clear_user_flag() const;
1197 
1203  void
1204  recursively_set_user_flag() const;
1205 
1211  void
1212  recursively_clear_user_flag() const;
1213 
1219  void
1220  clear_user_data() const;
1221 
1233  void
1234  set_user_pointer(void *p) const;
1235 
1241  void
1242  clear_user_pointer() const;
1243 
1259  void *
1260  user_pointer() const;
1261 
1283  void
1284  recursively_set_user_pointer(void *p) const;
1285 
1292  void
1293  recursively_clear_user_pointer() const;
1294 
1304  void
1305  set_user_index(const unsigned int p) const;
1306 
1312  void
1313  clear_user_index() const;
1314 
1326  unsigned int
1327  user_index() const;
1328 
1346  void
1347  recursively_set_user_index(const unsigned int p) const;
1348 
1357  void
1358  recursively_clear_user_index() const;
1377  double
1378  diameter() const;
1379 
1406  std::pair<Point<spacedim>, double>
1407  enclosing_ball() const;
1408 
1413  bounding_box() const;
1414 
1424  double
1425  extent_in_direction(const unsigned int axis) const;
1426 
1430  double
1431  minimum_vertex_distance() const;
1432 
1447  intermediate_point(const Point<structdim> &coordinates) const;
1448 
1472  real_to_unit_cell_affine_approximation(const Point<spacedim> &point) const;
1473 
1509  center(const bool respect_manifold = false,
1510  const bool interpolate_from_surrounding = false) const;
1511 
1530  barycenter() const;
1531 
1557  double
1558  measure() const;
1559 
1574  bool
1575  is_translation_of(
1577 
1583 private:
1588  void
1589  set_boundary_id_internal(const types::boundary_id id) const;
1590 
1595  void
1596  set(const ::internal::TriangulationImplementation::TriaObject<structdim>
1597  &o) const;
1598 
1606  void
1607  set_line_orientation(const unsigned int line, const bool orientation) const;
1608 
1619  void
1620  set_face_orientation(const unsigned int face, const bool orientation) const;
1621 
1628  void
1629  set_face_flip(const unsigned int face, const bool flip) const;
1630 
1637  void
1638  set_face_rotation(const unsigned int face, const bool rotation) const;
1639 
1643  void
1644  set_used_flag() const;
1645 
1649  void
1650  clear_used_flag() const;
1651 
1660  void
1661  set_refinement_case(const RefinementCase<structdim> &ref_case) const;
1662 
1670  void
1671  clear_refinement_case() const;
1672 
1679  void
1680  set_children(const unsigned int i, const int index) const;
1681 
1686  void
1687  clear_children() const;
1688 
1689 private:
1690  template <int, int>
1691  friend class Triangulation;
1692 
1693  friend struct ::internal::TriangulationImplementation::Implementation;
1694  friend struct ::internal::TriaAccessorImplementation::Implementation;
1695 };
1696 
1697 
1698 
1718 template <int dim, int spacedim>
1719 class TriaAccessor<0, dim, spacedim>
1720 {
1721 public:
1727  static const unsigned int space_dimension = spacedim;
1728 
1734  static const unsigned int dimension = dim;
1735 
1741  static const unsigned int structure_dimension = 0;
1742 
1746  using AccessorData = void;
1747 
1753  const unsigned int vertex_index);
1754 
1760  TriaAccessor(const Triangulation<dim, spacedim> *tria = nullptr,
1761  const int level = 0,
1762  const int index = 0,
1763  const AccessorData * = nullptr);
1764 
1768  template <int structdim2, int dim2, int spacedim2>
1770 
1774  template <int structdim2, int dim2, int spacedim2>
1776 
1781  state() const;
1782 
1787  static int
1788  level();
1789 
1794  int
1795  index() const;
1796 
1806  void
1807  operator++();
1808 
1812  void
1813  operator--();
1817  bool
1818  operator==(const TriaAccessor &) const;
1819 
1823  bool
1824  operator!=(const TriaAccessor &) const;
1825 
1853  unsigned int
1854  vertex_index(const unsigned int i = 0) const;
1855 
1861  Point<spacedim> &
1862  vertex(const unsigned int i = 0) const;
1863 
1868  typename ::internal::TriangulationImplementation::
1869  Iterators<dim, spacedim>::line_iterator static line(const unsigned int);
1870 
1874  static unsigned int
1875  line_index(const unsigned int i);
1876 
1880  static typename ::internal::TriangulationImplementation::
1881  Iterators<dim, spacedim>::quad_iterator
1882  quad(const unsigned int i);
1883 
1887  static unsigned int
1888  quad_index(const unsigned int i);
1889 
1905  double
1906  diameter() const;
1907 
1915  double
1916  extent_in_direction(const unsigned int axis) const;
1917 
1926  center(const bool respect_manifold = false,
1927  const bool interpolate_from_surrounding = false) const;
1928 
1937  double
1938  measure() const;
1953  static bool
1954  face_orientation(const unsigned int face);
1955 
1959  static bool
1960  face_flip(const unsigned int face);
1961 
1965  static bool
1966  face_rotation(const unsigned int face);
1967 
1971  static bool
1972  line_orientation(const unsigned int line);
1973 
1988  static bool
1989  has_children();
1990 
1995  static unsigned int
1996  n_children();
1997 
2002  static unsigned int
2003  number_of_children();
2004 
2008  static unsigned int
2009  max_refinement_depth();
2010 
2015  child(const unsigned int);
2016 
2021  isotropic_child(const unsigned int);
2022 
2026  static RefinementCase<0>
2027  refinement_case();
2028 
2032  static int
2033  child_index(const unsigned int i);
2034 
2038  static int
2039  isotropic_child_index(const unsigned int i);
2047  bool
2048  used() const;
2049 
2050 protected:
2058  void
2059  copy_from(const TriaAccessor &);
2060 
2065 
2069  unsigned int global_vertex_index;
2070 
2071 private:
2072  template <typename Accessor>
2073  friend class TriaRawIterator;
2074  template <typename Accessor>
2075  friend class TriaIterator;
2076  template <typename Accessor>
2077  friend class TriaActiveIterator;
2078 };
2079 
2080 
2081 
2099 template <int spacedim>
2100 class TriaAccessor<0, 1, spacedim>
2101 {
2102 public:
2108  static const unsigned int space_dimension = spacedim;
2109 
2115  static const unsigned int dimension = 1;
2116 
2122  static const unsigned int structure_dimension = 0;
2123 
2127  using AccessorData = void;
2128 
2134  {
2146  right_vertex
2147  };
2148 
2161  const VertexKind vertex_kind,
2162  const unsigned int vertex_index);
2163 
2169  TriaAccessor(const Triangulation<1, spacedim> *tria = nullptr,
2170  const int = 0,
2171  const int = 0,
2172  const AccessorData * = nullptr);
2173 
2177  template <int structdim2, int dim2, int spacedim2>
2179 
2183  template <int structdim2, int dim2, int spacedim2>
2185 
2190  void
2191  copy_from(const TriaAccessor &);
2192 
2199  state();
2200 
2205  static int
2206  level();
2207 
2212  int
2213  index() const;
2214 
2225  void
2226  operator++() const;
2227 
2232  void
2233  operator--() const;
2237  bool
2238  operator==(const TriaAccessor &) const;
2239 
2243  bool
2244  operator!=(const TriaAccessor &) const;
2245 
2272  unsigned int
2273  vertex_index(const unsigned int i = 0) const;
2274 
2280  Point<spacedim> &
2281  vertex(const unsigned int i = 0) const;
2282 
2288  center() const;
2289 
2294  typename ::internal::TriangulationImplementation::
2295  Iterators<1, spacedim>::line_iterator static line(const unsigned int);
2296 
2303  static unsigned int
2304  line_index(const unsigned int i);
2305 
2309  static typename ::internal::TriangulationImplementation::
2310  Iterators<1, spacedim>::quad_iterator
2311  quad(const unsigned int i);
2312 
2319  static unsigned int
2320  quad_index(const unsigned int i);
2321 
2331  bool
2332  at_boundary() const;
2333 
2349  boundary_id() const;
2350 
2354  const Manifold<1, spacedim> &
2355  get_manifold() const;
2356 
2364  manifold_id() const;
2365 
2366 
2377  static bool
2378  face_orientation(const unsigned int face);
2379 
2383  static bool
2384  face_flip(const unsigned int face);
2385 
2389  static bool
2390  face_rotation(const unsigned int face);
2391 
2395  static bool
2396  line_orientation(const unsigned int line);
2397 
2412  static bool
2413  has_children();
2414 
2419  static unsigned int
2420  n_children();
2421 
2426  static unsigned int
2427  number_of_children();
2428 
2432  static unsigned int
2433  max_refinement_depth();
2434 
2439  child(const unsigned int);
2440 
2445  isotropic_child(const unsigned int);
2446 
2450  static RefinementCase<0>
2451  refinement_case();
2452 
2456  static int
2457  child_index(const unsigned int i);
2458 
2462  static int
2463  isotropic_child_index(const unsigned int i);
2496  void
2497  set_boundary_id(const types::boundary_id);
2498 
2505  void
2506  set_manifold_id(const types::manifold_id);
2507 
2519  void
2520  set_all_boundary_ids(const types::boundary_id);
2521 
2533  void
2534  set_all_manifold_ids(const types::manifold_id);
2542  bool
2543  used() const;
2544 
2545 protected:
2550 
2556 
2560  unsigned int global_vertex_index;
2561 };
2562 
2563 
2564 
2581 template <int dim, int spacedim = dim>
2582 class CellAccessor : public TriaAccessor<dim, dim, spacedim>
2583 {
2584 public:
2589 
2594 
2605  CellAccessor(const Triangulation<dim, spacedim> *parent = nullptr,
2606  const int level = -1,
2607  const int index = -1,
2608  const AccessorData * local_data = nullptr);
2609 
2613  CellAccessor(const TriaAccessor<dim, dim, spacedim> &cell_accessor);
2614 
2627  template <int structdim2, int dim2, int spacedim2>
2629 
2634  template <int structdim2, int dim2, int spacedim2>
2636 
2646  void
2647  operator=(const CellAccessor<dim, spacedim> &) = delete;
2648 
2665  child(const unsigned int i) const;
2666 
2670  TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>
2671  face(const unsigned int i) const;
2672 
2682  unsigned int
2683  face_index(const unsigned int i) const;
2684 
2733  neighbor_child_on_subface(const unsigned int face_no,
2734  const unsigned int subface_no) const;
2735 
2758  neighbor(const unsigned int i) const;
2759 
2764  int
2765  neighbor_index(const unsigned int i) const;
2766 
2771  int
2772  neighbor_level(const unsigned int i) const;
2773 
2785  unsigned int
2786  neighbor_of_neighbor(const unsigned int neighbor) const;
2787 
2798  bool
2799  neighbor_is_coarser(const unsigned int neighbor) const;
2800 
2815  std::pair<unsigned int, unsigned int>
2816  neighbor_of_coarser_neighbor(const unsigned int neighbor) const;
2817 
2824  unsigned int
2825  neighbor_face_no(const unsigned int neighbor) const;
2826 
2830  static bool
2831  is_level_cell();
2832 
2846  bool
2847  has_periodic_neighbor(const unsigned int i) const;
2848 
2866  periodic_neighbor(const unsigned int i) const;
2867 
2876  neighbor_or_periodic_neighbor(const unsigned int i) const;
2877 
2893  periodic_neighbor_child_on_subface(const unsigned int face_no,
2894  const unsigned int subface_no) const;
2895 
2906  std::pair<unsigned int, unsigned int>
2907  periodic_neighbor_of_coarser_periodic_neighbor(const unsigned face_no) const;
2908 
2914  int
2915  periodic_neighbor_index(const unsigned int i) const;
2916 
2922  int
2923  periodic_neighbor_level(const unsigned int i) const;
2924 
2939  unsigned int
2940  periodic_neighbor_of_periodic_neighbor(const unsigned int i) const;
2941 
2947  unsigned int
2948  periodic_neighbor_face_no(const unsigned int i) const;
2949 
2956  bool
2957  periodic_neighbor_is_coarser(const unsigned int i) const;
2958 
2975  bool
2976  at_boundary(const unsigned int i) const;
2977 
2986  bool
2987  at_boundary() const;
2988 
2996  bool
2997  has_boundary_lines() const;
3024  refine_flag_set() const;
3025 
3043  void
3044  set_refine_flag(const RefinementCase<dim> ref_case =
3046 
3050  void
3051  clear_refine_flag() const;
3052 
3060  bool
3061  flag_for_face_refinement(
3062  const unsigned int face_no,
3063  const RefinementCase<dim - 1> &face_refinement_case =
3065 
3071  bool
3072  flag_for_line_refinement(const unsigned int line_no) const;
3073 
3083  subface_case(const unsigned int face_no) const;
3084 
3088  bool
3089  coarsen_flag_set() const;
3090 
3095  void
3096  set_coarsen_flag() const;
3097 
3101  void
3102  clear_coarsen_flag() const;
3126  material_id() const;
3127 
3139  void
3140  set_material_id(const types::material_id new_material_id) const;
3141 
3150  void
3151  recursively_set_material_id(const types::material_id new_material_id) const;
3178  subdomain_id() const;
3179 
3195  void
3196  set_subdomain_id(const types::subdomain_id new_subdomain_id) const;
3197 
3203  level_subdomain_id() const;
3204 
3209  void
3210  set_level_subdomain_id(
3211  const types::subdomain_id new_level_subdomain_id) const;
3212 
3213 
3229  void
3230  recursively_set_subdomain_id(
3231  const types::subdomain_id new_subdomain_id) const;
3249  bool
3250  direction_flag() const;
3251 
3277  unsigned int
3278  active_cell_index() const;
3279 
3287  int
3288  parent_index() const;
3289 
3296  parent() const;
3297 
3317  bool
3318  active() const;
3319 
3339  bool
3340  is_locally_owned() const;
3341 
3346  bool
3347  is_locally_owned_on_level() const;
3348 
3372  bool
3373  is_ghost() const;
3374 
3401  bool
3402  is_artificial() const;
3403 
3417  bool
3418  point_inside(const Point<spacedim> &p) const;
3419 
3428  void
3429  set_neighbor(const unsigned int i,
3430  const TriaIterator<CellAccessor<dim, spacedim>> &pointer) const;
3431 
3445  CellId
3446  id() const;
3447 
3456  DeclException0(ExcRefineCellNotActive);
3460  DeclException0(ExcCellFlaggedForRefinement);
3464  DeclException0(ExcCellFlaggedForCoarsening);
3465 
3466 protected:
3482  unsigned int
3483  neighbor_of_neighbor_internal(const unsigned int neighbor) const;
3484 
3490  template <int dim_, int spacedim_>
3491  bool
3492  point_inside_codim(const Point<spacedim_> &p) const;
3493 
3494 
3495 
3496 private:
3501  void
3502  set_active_cell_index(const unsigned int active_cell_index);
3503 
3507  void
3508  set_parent(const unsigned int parent_index);
3509 
3516  void
3517  set_direction_flag(const bool new_direction_flag) const;
3518 
3519  template <int, int>
3520  friend class Triangulation;
3521 
3522  friend struct ::internal::TriangulationImplementation::Implementation;
3523 };
3524 
3525 
3526 
3527 /* ----- declaration of explicit specializations and general templates ----- */
3528 
3529 
3530 template <int structdim, int dim, int spacedim>
3531 template <typename OtherAccessor>
3533  const OtherAccessor &)
3534 {
3535  Assert(false,
3536  ExcMessage("You are attempting an illegal conversion between "
3537  "iterator/accessor types. The constructor you call "
3538  "only exists to make certain template constructs "
3539  "easier to write as dimension independent code but "
3540  "the conversion is not valid in the current context."));
3541 }
3542 
3543 
3544 
3545 template <int structdim, int dim, int spacedim>
3546 template <int structdim2, int dim2, int spacedim2>
3549 {
3550  Assert(false,
3551  ExcMessage("You are attempting an illegal conversion between "
3552  "iterator/accessor types. The constructor you call "
3553  "only exists to make certain template constructs "
3554  "easier to write as dimension independent code but "
3555  "the conversion is not valid in the current context."));
3556 }
3557 
3558 
3559 
3560 template <int dim, int spacedim>
3561 template <int structdim2, int dim2, int spacedim2>
3564 {
3565  Assert(false,
3566  ExcMessage("You are attempting an illegal conversion between "
3567  "iterator/accessor types. The constructor you call "
3568  "only exists to make certain template constructs "
3569  "easier to write as dimension independent code but "
3570  "the conversion is not valid in the current context."));
3571 }
3572 
3573 
3574 
3575 template <int structdim, int dim, int spacedim>
3576 template <int structdim2, int dim2, int spacedim2>
3579 {
3580  Assert(false,
3581  ExcMessage("You are attempting an illegal conversion between "
3582  "iterator/accessor types. The constructor you call "
3583  "only exists to make certain template constructs "
3584  "easier to write as dimension independent code but "
3585  "the conversion is not valid in the current context."));
3586 }
3587 
3588 
3589 
3590 template <int dim, int spacedim>
3591 template <int structdim2, int dim2, int spacedim2>
3594 {
3595  Assert(false,
3596  ExcMessage("You are attempting an illegal conversion between "
3597  "iterator/accessor types. The constructor you call "
3598  "only exists to make certain template constructs "
3599  "easier to write as dimension independent code but "
3600  "the conversion is not valid in the current context."));
3601 }
3602 
3603 
3604 #ifndef DOXYGEN
3605 
3606 template <>
3607 bool
3609 template <>
3610 bool
3612 template <>
3613 bool
3615 template <>
3616 bool
3618 template <>
3619 bool
3621 template <>
3622 bool
3624 // -------------------------------------------------------------------
3625 
3626 template <>
3627 void
3629 
3630 #endif // DOXYGEN
3631 
3632 DEAL_II_NAMESPACE_CLOSE
3633 
3634 // include more templates in debug and optimized mode
3635 #include "tria_accessor.templates.h"
3636 
3637 #endif
unsigned int manifold_id
Definition: types.h:123
static::ExceptionBase & ExcCellHasNoChildren()
InvalidAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
static::ExceptionBase & ExcNeighborIsNotCoarser()
const Triangulation< dim, spacedim > * tria
static::ExceptionBase & ExcNeighborIsCoarser()
static::ExceptionBase & ExcDereferenceInvalidObject(AccessorType arg1)
const Triangulation< 1, spacedim > * tria
unsigned int material_id
Definition: types.h:134
static::ExceptionBase & ExcFacesHaveNoLevel()
static::ExceptionBase & ExcCellNotUsed()
static::ExceptionBase & ExcMessage(std::string arg1)
unsigned int subdomain_id
Definition: types.h:43
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:408
const Triangulation< dim, spacedim > * tria
static::ExceptionBase & ExcSetOnlyEvenChildren(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1227
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:397
#define DeclException0(Exception0)
Definition: exceptions.h:385
static::ExceptionBase & ExcCellNotActive()
static::ExceptionBase & ExcNoPeriodicNeighbor()
static::ExceptionBase & ExcCantSetChildren(int arg1)
bool point_inside(const Point< spacedim > &p) const
Definition: cell_id.h:63
CellAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
static::ExceptionBase & ExcCellHasNoParent()
TriaAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
Iterator reached end of container.
Iterator points to a valid object.
void set_all_manifold_ids(const types::manifold_id) const
unsigned int boundary_id
Definition: types.h:111
typename::internal::TriaAccessorImplementation::PresentLevelType< structdim, dim >::type present_level
static::ExceptionBase & ExcInternalError()
static::ExceptionBase & ExcCantCompareIterators()