Reference documentation for deal.II version 9.3.2
\(\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\}}\)
tria_accessor.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2021 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 
25 #include <deal.II/base/point.h>
26 
27 #include <deal.II/grid/cell_id.h>
31 
33 #include <boost/container/small_vector.hpp>
35 
36 #include <utility>
37 
38 
40 
41 // Forward declarations
42 #ifndef DOXYGEN
43 template <int dim, int spacedim>
44 class Triangulation;
45 template <typename Accessor>
46 class TriaRawIterator;
47 template <typename Accessor>
48 class TriaIterator;
49 template <typename Accessor>
50 class TriaActiveIterator;
51 
52 namespace parallel
53 {
54  template <int dim, int spacedim>
55  class TriangulationBase;
56 }
57 
58 template <int dim, int spacedim>
59 class Manifold;
60 
61 template <int dim, int spacedim>
62 class Mapping;
63 #endif
64 
65 namespace internal
66 {
67  namespace TriangulationImplementation
68  {
69  class TriaObjects;
70  struct Implementation;
71  struct ImplementationMixedMesh;
72  } // namespace TriangulationImplementation
73 
74  namespace TriaAccessorImplementation
75  {
76  struct Implementation;
77 
83  template <int structdim, int dim>
85  {
86  struct type
87  {
91  type() = default;
92 
96  type(const int level)
97  {
98  Assert(level == 0, ExcInternalError());
99  (void)level; // removes -Wunused-parameter warning in optimized mode
100  }
101 
105  operator int() const
106  {
107  return 0;
108  }
109 
110  void
111  operator++() const
112  {
113  Assert(false, ExcInternalError());
114  }
115 
116  void
117  operator--() const
118  {
119  Assert(false, ExcInternalError());
120  }
121  };
122  };
123 
124 
130  template <int dim>
131  struct PresentLevelType<dim, dim>
132  {
133  using type = int;
134  };
135  } // namespace TriaAccessorImplementation
136 } // namespace internal
137 template <int structdim, int dim, int spacedim>
138 class TriaAccessor;
139 template <int dim, int spacedim>
140 class TriaAccessor<0, dim, spacedim>;
141 template <int spacedim>
142 class TriaAccessor<0, 1, spacedim>;
143 
144 // note: the file tria_accessor.templates.h is included at the end of
145 // this file. this includes a lot of templates. originally, this was
146 // only done in debug mode, but led to cyclic reduction problems and
147 // so is now on by default.
148 
149 
154 {
159  "The operation you are attempting can only be performed for "
160  "(cell, face, or edge) iterators that point to valid "
161  "objects. These objects need not necessarily be active, "
162  "i.e., have no children, but they need to be part of a "
163  "triangulation. (The objects pointed to by an iterator "
164  "may -- after coarsening -- also be objects that used "
165  "to be part of a triangulation, but are now no longer "
166  "used. Their memory location may have been retained "
167  "for re-use upon the next mesh refinement, but is "
168  "currently unused.)");
179  "The operation you are attempting can only be performed for "
180  "(cell, face, or edge) iterators that point to 'active' "
181  "objects. 'Active' objects are those that do not have "
182  "children (in the case of cells), or that are part of "
183  "an active cell (in the case of faces or edges). However, "
184  "the object on which you are trying the current "
185  "operation is not 'active' in this sense.");
192  "The operation you are attempting can only be performed for "
193  "(cell, face, or edge) iterators that have children, "
194  "but the object on which you are trying the current "
195  "operation does not have any.");
203  "The operation you are attempting can only be performed for "
204  "(cell, face, or edge) iterators that have a parent object, "
205  "but the object on which you are trying the current "
206  "operation does not have one -- i.e., it is on the "
207  "coarsest level of the triangulation.");
212  int,
213  << "You can only set the child index if the cell does not "
214  << "currently have children registered; or you can clear it. "
215  << "The given index was " << arg1
216  << " (-1 means: clear children).");
220  template <typename AccessorType>
222  AccessorType,
223  << "You tried to dereference an iterator for which this "
224  << "is not possible. More information on this iterator: "
225  << "index=" << arg1.index() << ", state="
226  << (arg1.state() == IteratorState::valid ?
227  "valid" :
228  (arg1.state() == IteratorState::past_the_end ?
229  "past_the_end" :
230  "invalid")));
235  "Iterators can only be compared if they point to the same "
236  "triangulation, or if neither of them are associated "
237  "with a triangulation.");
238  // TODO: Write documentation!
243  // TODO: Write documentation!
263  // TODO: Write documentation!
269  int,
270  << "You can only set the child index of an even numbered child."
271  << "The number of the child given was " << arg1 << ".");
272 } // namespace TriaAccessorExceptions
273 
274 
300 template <int structdim, int dim, int spacedim = dim>
302 {
303 public:
309  static const unsigned int space_dimension = spacedim;
310 
316  static const unsigned int dimension = dim;
317 
323  static const unsigned int structure_dimension = structdim;
324 
334  void
335  operator=(const TriaAccessorBase *) = delete;
336 
337 protected:
343  using AccessorData = void;
344 
349  const int level = -1,
350  const int index = -1,
351  const AccessorData * = nullptr);
352 
357 
365  void
367 
373 
384  bool
385  operator<(const TriaAccessorBase &other) const;
386 
387 protected:
391  bool
392  operator==(const TriaAccessorBase &) const;
393 
397  bool
398  operator!=(const TriaAccessorBase &) const;
399 
413  void
415 
423  void
433  objects() const;
434 
435 public:
441  using LocalData = void *;
442 
466  int
467  level() const;
468 
495  int
496  index() const;
497 
503  state() const;
504 
511 
515 protected:
520  typename ::internal::TriaAccessorImplementation::
521  PresentLevelType<structdim, dim>::type present_level;
522 
528 
533 
534 private:
535  template <typename Accessor>
536  friend class TriaRawIterator;
537  template <typename Accessor>
538  friend class TriaIterator;
539  template <typename Accessor>
540  friend class TriaActiveIterator;
541 };
542 
543 
544 
565 template <int structdim, int dim, int spacedim = dim>
566 class InvalidAccessor : public TriaAccessorBase<structdim, dim, spacedim>
567 {
568 public:
572  using AccessorData =
574 
583  const int level = -1,
584  const int index = -1,
585  const AccessorData * local_data = nullptr);
586 
595 
600  template <typename OtherAccessor>
601  InvalidAccessor(const OtherAccessor &);
602 
606  void
608 
612  bool
613  operator==(const InvalidAccessor &) const;
614  bool
615  operator!=(const InvalidAccessor &) const;
616 
620  void
621  operator++() const;
622  void
623  operator--() const;
624 
629  bool
630  used() const;
631 
636  bool
637  has_children() const;
638 
643  manifold_id() const;
644 
648  unsigned int
649  user_index() const;
650 
654  void
655  set_user_index(const unsigned int p) const;
656 
660  void
662 
667  vertex(const unsigned int i) const;
668 
673  typename ::internal::TriangulationImplementation::
674  Iterators<dim, spacedim>::line_iterator
675  line(const unsigned int i) const;
676 
681  typename ::internal::TriangulationImplementation::
682  Iterators<dim, spacedim>::quad_iterator
683  quad(const unsigned int i) const;
684 };
685 
686 
687 
705 template <int structdim, int dim, int spacedim>
706 class TriaAccessor : public TriaAccessorBase<structdim, dim, spacedim>
707 {
708 public:
712  using AccessorData =
714 
719  const int level = -1,
720  const int index = -1,
721  const AccessorData * local_data = nullptr);
722 
727  TriaAccessor(const TriaAccessor &) = default;
728 
732  TriaAccessor(TriaAccessor &&) = default; // NOLINT
733 
746  template <int structdim2, int dim2, int spacedim2>
748 
753  template <int structdim2, int dim2, int spacedim2>
755 
765  TriaAccessor &
766  operator=(const TriaAccessor &) = delete;
767 
771  TriaAccessor &
772  operator=(TriaAccessor &&) = default; // NOLINT
773 
777  ~TriaAccessor() = default;
778 
785  bool
786  used() const;
787 
800  vertex_iterator(const unsigned int i) const;
801 
817  unsigned int
818  vertex_index(const unsigned int i) const;
819 
858  vertex(const unsigned int i) const;
859 
863  typename ::internal::TriangulationImplementation::
864  Iterators<dim, spacedim>::line_iterator
865  line(const unsigned int i) const;
866 
873  unsigned int
874  line_index(const unsigned int i) const;
875 
879  typename ::internal::TriangulationImplementation::
880  Iterators<dim, spacedim>::quad_iterator
881  quad(const unsigned int i) const;
882 
889  unsigned int
890  quad_index(const unsigned int i) const;
913  bool
914  face_orientation(const unsigned int face) const;
915 
925  bool
926  face_flip(const unsigned int face) const;
927 
937  bool
938  face_rotation(const unsigned int face) const;
939 
952  bool
953  line_orientation(const unsigned int line) const;
968  bool
969  has_children() const;
970 
975  unsigned int
976  n_children() const;
977 
982  unsigned int
984 
998  unsigned int
1000 
1014  unsigned int
1016 
1021  child(const unsigned int i) const;
1022 
1027  unsigned int
1030 
1040  isotropic_child(const unsigned int i) const;
1041 
1047 
1053  int
1054  child_index(const unsigned int i) const;
1055 
1061  int
1062  isotropic_child_index(const unsigned int i) const;
1085  boundary_id() const;
1086 
1116  void
1118 
1149  void
1151 
1159  bool
1160  at_boundary() const;
1161 
1171  const Manifold<dim, spacedim> &
1172  get_manifold() const;
1173 
1195  manifold_id() const;
1196 
1214  void
1216 
1230  void
1232 
1249  bool
1250  user_flag_set() const;
1251 
1257  void
1258  set_user_flag() const;
1259 
1265  void
1267 
1273  void
1275 
1281  void
1283 
1289  void
1291 
1303  void
1304  set_user_pointer(void *p) const;
1305 
1311  void
1313 
1329  void *
1330  user_pointer() const;
1331 
1353  void
1355 
1362  void
1364 
1374  void
1375  set_user_index(const unsigned int p) const;
1376 
1382  void
1384 
1396  unsigned int
1397  user_index() const;
1398 
1416  void
1417  recursively_set_user_index(const unsigned int p) const;
1418 
1427  void
1463  double
1464  diameter() const;
1465 
1492  std::pair<Point<spacedim>, double>
1494 
1504  bounding_box() const;
1505 
1515  double
1516  extent_in_direction(const unsigned int axis) const;
1517 
1521  double
1523 
1538  intermediate_point(const Point<structdim> &coordinates) const;
1539 
1564 
1600  center(const bool respect_manifold = false,
1601  const bool interpolate_from_surrounding = false) const;
1602 
1621  barycenter() const;
1622 
1648  double
1649  measure() const;
1650 
1665  bool
1668 
1674 
1678  unsigned int
1679  n_vertices() const;
1680 
1684  unsigned int
1685  n_lines() const;
1686 
1692  unsigned int
1693  n_faces() const;
1694 
1701 
1707  line_indices() const;
1708 
1716  face_indices() const;
1717 
1722 private:
1727  void
1729 
1737  void
1739  const std::initializer_list<int> &new_indices) const;
1740 
1744  void
1746  const std::initializer_list<unsigned int> &new_indices) const;
1747 
1755  void
1756  set_line_orientation(const unsigned int line, const bool orientation) const;
1757 
1768  void
1769  set_face_orientation(const unsigned int face, const bool orientation) const;
1770 
1777  void
1778  set_face_flip(const unsigned int face, const bool flip) const;
1779 
1786  void
1787  set_face_rotation(const unsigned int face, const bool rotation) const;
1788 
1792  void
1793  set_used_flag() const;
1794 
1798  void
1800 
1809  void
1811 
1819  void
1821 
1828  void
1829  set_children(const unsigned int i, const int index) const;
1830 
1835  void
1837 
1838 private:
1839  template <int, int>
1840  friend class Triangulation;
1841 
1842  friend struct ::internal::TriangulationImplementation::Implementation;
1843  friend struct ::internal::TriangulationImplementation::
1844  ImplementationMixedMesh;
1845  friend struct ::internal::TriaAccessorImplementation::Implementation;
1846 };
1847 
1848 
1849 
1868 template <int dim, int spacedim>
1869 class TriaAccessor<0, dim, spacedim>
1870 {
1871 public:
1877  static const unsigned int space_dimension = spacedim;
1878 
1884  static const unsigned int dimension = dim;
1885 
1891  static const unsigned int structure_dimension = 0;
1892 
1896  using AccessorData = void;
1897 
1903  const unsigned int vertex_index);
1904 
1911  const int level = 0,
1912  const int index = 0,
1913  const AccessorData * = nullptr);
1914 
1918  template <int structdim2, int dim2, int spacedim2>
1920 
1924  template <int structdim2, int dim2, int spacedim2>
1926 
1931  state() const;
1932 
1937  static int
1939 
1944  int
1945  index() const;
1946 
1953 
1963  void
1965 
1969  void
1974  bool
1975  operator==(const TriaAccessor &) const;
1976 
1980  bool
1981  operator!=(const TriaAccessor &) const;
1982 
2010  unsigned int
2011  vertex_index(const unsigned int i = 0) const;
2012 
2018  Point<spacedim> &
2019  vertex(const unsigned int i = 0) const;
2020 
2025  typename ::internal::TriangulationImplementation::
2026  Iterators<dim, spacedim>::line_iterator static line(const unsigned int);
2027 
2031  static unsigned int
2032  line_index(const unsigned int i);
2033 
2037  static typename ::internal::TriangulationImplementation::
2038  Iterators<dim, spacedim>::quad_iterator
2039  quad(const unsigned int i);
2040 
2044  static unsigned int
2045  quad_index(const unsigned int i);
2046 
2062  double
2063  diameter() const;
2064 
2072  double
2073  extent_in_direction(const unsigned int axis) const;
2074 
2083  center(const bool respect_manifold = false,
2084  const bool interpolate_from_surrounding = false) const;
2085 
2094  double
2095  measure() const;
2110  static bool
2111  face_orientation(const unsigned int face);
2112 
2116  static bool
2117  face_flip(const unsigned int face);
2118 
2122  static bool
2123  face_rotation(const unsigned int face);
2124 
2128  static bool
2129  line_orientation(const unsigned int line);
2130 
2145  static bool
2147 
2152  static unsigned int
2154 
2159  static unsigned int
2161 
2166  static unsigned int
2168 
2169 
2173  static unsigned int
2175 
2179  static unsigned int
2181 
2186  child(const unsigned int);
2187 
2192  isotropic_child(const unsigned int);
2193 
2197  static RefinementCase<0>
2199 
2203  static int
2204  child_index(const unsigned int i);
2205 
2209  static int
2210  isotropic_child_index(const unsigned int i);
2218  bool
2219  used() const;
2220 
2221 protected:
2229  void
2231 
2240  bool
2241  operator<(const TriaAccessor &other) const;
2242 
2247 
2251  unsigned int global_vertex_index;
2252 
2253 private:
2254  template <typename Accessor>
2255  friend class TriaRawIterator;
2256  template <typename Accessor>
2257  friend class TriaIterator;
2258  template <typename Accessor>
2259  friend class TriaActiveIterator;
2260 };
2261 
2262 
2263 
2280 template <int spacedim>
2281 class TriaAccessor<0, 1, spacedim>
2282 {
2283 public:
2289  static const unsigned int space_dimension = spacedim;
2290 
2296  static const unsigned int dimension = 1;
2297 
2303  static const unsigned int structure_dimension = 0;
2304 
2308  using AccessorData = void;
2309 
2315  {
2327  right_vertex
2328  };
2329 
2342  const VertexKind vertex_kind,
2343  const unsigned int vertex_index);
2344 
2351  const int = 0,
2352  const int = 0,
2353  const AccessorData * = nullptr);
2354 
2358  template <int structdim2, int dim2, int spacedim2>
2360 
2364  template <int structdim2, int dim2, int spacedim2>
2366 
2371  void
2373 
2381 
2386  static int
2388 
2393  int
2394  index() const;
2395 
2402 
2413  void
2414  operator++() const;
2415 
2420  void
2421  operator--() const;
2425  bool
2426  operator==(const TriaAccessor &) const;
2427 
2431  bool
2432  operator!=(const TriaAccessor &) const;
2433 
2442  bool
2443  operator<(const TriaAccessor &other) const;
2444 
2471  unsigned int
2472  vertex_index(const unsigned int i = 0) const;
2473 
2479  Point<spacedim> &
2480  vertex(const unsigned int i = 0) const;
2481 
2487  center() const;
2488 
2493  typename ::internal::TriangulationImplementation::
2494  Iterators<1, spacedim>::line_iterator static line(const unsigned int);
2495 
2502  static unsigned int
2503  line_index(const unsigned int i);
2504 
2508  static typename ::internal::TriangulationImplementation::
2509  Iterators<1, spacedim>::quad_iterator
2510  quad(const unsigned int i);
2511 
2518  static unsigned int
2519  quad_index(const unsigned int i);
2520 
2530  bool
2531  at_boundary() const;
2532 
2548  boundary_id() const;
2549 
2553  const Manifold<1, spacedim> &
2554  get_manifold() const;
2555 
2563  manifold_id() const;
2564 
2565 
2576  static bool
2577  face_orientation(const unsigned int face);
2578 
2582  static bool
2583  face_flip(const unsigned int face);
2584 
2588  static bool
2589  face_rotation(const unsigned int face);
2590 
2594  static bool
2595  line_orientation(const unsigned int line);
2596 
2611  static bool
2613 
2618  static unsigned int
2620 
2625  static unsigned int
2627 
2632  static unsigned int
2634 
2635 
2639  static unsigned int
2641 
2645  static unsigned int
2647 
2652  child(const unsigned int);
2653 
2658  isotropic_child(const unsigned int);
2659 
2663  static RefinementCase<0>
2665 
2669  static int
2670  child_index(const unsigned int i);
2671 
2675  static int
2676  isotropic_child_index(const unsigned int i);
2709  void
2711 
2718  void
2720 
2732  void
2734 
2746  void
2755  bool
2756  used() const;
2757 
2763 
2767  unsigned int
2768  n_vertices() const;
2769 
2773  unsigned int
2774  n_lines() const;
2775 
2782 
2788  line_indices() const;
2789 
2790 protected:
2795 
2801 
2805  unsigned int global_vertex_index;
2806 };
2807 
2808 
2809 
2825 template <int dim, int spacedim = dim>
2826 class CellAccessor : public TriaAccessor<dim, dim, spacedim>
2827 {
2828 public:
2833 
2838 
2850  const int level = -1,
2851  const int index = -1,
2852  const AccessorData * local_data = nullptr);
2853 
2858 
2871  template <int structdim2, int dim2, int spacedim2>
2873 
2878  template <int structdim2, int dim2, int spacedim2>
2880 
2885 
2889  // NOLINTNEXTLINE OSX does not compile with noexcept
2891 
2895  ~CellAccessor() = default;
2896 
2908 
2912  // NOLINTNEXTLINE OSX does not compile with noexcept
2914  operator=(CellAccessor<dim, spacedim> &&) = default; // NOLINT
2915 
2932  child(const unsigned int i) const;
2933 
2937  boost::container::small_vector<TriaIterator<CellAccessor<dim, spacedim>>,
2940 
2944  TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>
2945  face(const unsigned int i) const;
2946 
2951  unsigned int
2954 
2958  boost::container::small_vector<
2959  TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>,
2962 
2972  unsigned int
2973  face_index(const unsigned int i) const;
2974 
3023  neighbor_child_on_subface(const unsigned int face_no,
3024  const unsigned int subface_no) const;
3025 
3051  neighbor(const unsigned int face_no) const;
3052 
3060  int
3061  neighbor_index(const unsigned int face_no) const;
3062 
3070  int
3071  neighbor_level(const unsigned int face_no) const;
3072 
3084  unsigned int
3085  neighbor_of_neighbor(const unsigned int face_no) const;
3086 
3097  bool
3098  neighbor_is_coarser(const unsigned int face_no) const;
3099 
3114  std::pair<unsigned int, unsigned int>
3115  neighbor_of_coarser_neighbor(const unsigned int neighbor) const;
3116 
3123  unsigned int
3124  neighbor_face_no(const unsigned int neighbor) const;
3125 
3129  static bool
3131 
3145  bool
3146  has_periodic_neighbor(const unsigned int i) const;
3147 
3165  periodic_neighbor(const unsigned int i) const;
3166 
3175  neighbor_or_periodic_neighbor(const unsigned int i) const;
3176 
3192  periodic_neighbor_child_on_subface(const unsigned int face_no,
3193  const unsigned int subface_no) const;
3194 
3205  std::pair<unsigned int, unsigned int>
3206  periodic_neighbor_of_coarser_periodic_neighbor(const unsigned face_no) const;
3207 
3213  int
3214  periodic_neighbor_index(const unsigned int i) const;
3215 
3221  int
3222  periodic_neighbor_level(const unsigned int i) const;
3223 
3238  unsigned int
3239  periodic_neighbor_of_periodic_neighbor(const unsigned int i) const;
3240 
3246  unsigned int
3247  periodic_neighbor_face_no(const unsigned int i) const;
3248 
3255  bool
3256  periodic_neighbor_is_coarser(const unsigned int i) const;
3257 
3274  bool
3275  at_boundary(const unsigned int i) const;
3276 
3285  bool
3286  at_boundary() const;
3287 
3295  bool
3296  has_boundary_lines() const;
3324 
3342  void
3345 
3349  void
3351 
3359  bool
3361  const unsigned int face_no,
3362  const RefinementCase<dim - 1> &face_refinement_case =
3364 
3370  bool
3371  flag_for_line_refinement(const unsigned int line_no) const;
3372 
3382  subface_case(const unsigned int face_no) const;
3383 
3387  bool
3389 
3394  void
3396 
3400  void
3425  material_id() const;
3426 
3438  void
3439  set_material_id(const types::material_id new_material_id) const;
3440 
3449  void
3450  recursively_set_material_id(const types::material_id new_material_id) const;
3477  subdomain_id() const;
3478 
3494  void
3495  set_subdomain_id(const types::subdomain_id new_subdomain_id) const;
3496 
3506  level_subdomain_id() const;
3507 
3512  void
3514  const types::subdomain_id new_level_subdomain_id) const;
3515 
3516 
3532  void
3534  const types::subdomain_id new_subdomain_id) const;
3552  global_active_cell_index() const;
3553 
3560  global_level_cell_index() const;
3561 
3575  bool
3576  direction_flag() const;
3577 
3603  unsigned int
3604  active_cell_index() const;
3605 
3613  int
3614  parent_index() const;
3615 
3622  parent() const;
3623 
3649  bool
3650  active() const;
3651 
3660  bool
3661  is_active() const;
3662 
3682  bool
3684 
3689  bool
3691 
3715  bool
3716  is_ghost() const;
3717 
3744  bool
3745  is_artificial() const;
3746 
3760  bool
3762 
3771  void
3772  set_neighbor(const unsigned int i,
3773  const TriaIterator<CellAccessor<dim, spacedim>> &pointer) const;
3774 
3788  CellId
3789  id() const;
3790 
3792 
3796  double
3797  diameter(const Mapping<dim, spacedim> &mapping) const;
3798 
3816 
3817 protected:
3833  unsigned int
3834  neighbor_of_neighbor_internal(const unsigned int neighbor) const;
3835 
3841  template <int dim_, int spacedim_>
3842  bool
3843  point_inside_codim(const Point<spacedim_> &p) const;
3844 
3845 
3846 
3847 private:
3852  void
3853  set_active_cell_index(const unsigned int active_cell_index) const;
3854 
3858  void
3860 
3864  void
3866 
3870  void
3871  set_parent(const unsigned int parent_index);
3872 
3879  void
3880  set_direction_flag(const bool new_direction_flag) const;
3881 
3882  template <int, int>
3883  friend class Triangulation;
3884 
3885  template <int, int>
3887 
3888  friend struct ::internal::TriangulationImplementation::Implementation;
3889  friend struct ::internal::TriangulationImplementation::
3890  ImplementationMixedMesh;
3891 };
3892 
3893 
3894 
3895 /* ----- declaration of explicit specializations and general templates ----- */
3896 
3897 
3898 template <int structdim, int dim, int spacedim>
3899 template <typename OtherAccessor>
3901  const OtherAccessor &)
3902 {
3903  Assert(false,
3904  ExcMessage("You are attempting an illegal conversion between "
3905  "iterator/accessor types. The constructor you call "
3906  "only exists to make certain template constructs "
3907  "easier to write as dimension independent code but "
3908  "the conversion is not valid in the current context."));
3909 }
3910 
3911 
3912 
3913 template <int structdim, int dim, int spacedim>
3914 template <int structdim2, int dim2, int spacedim2>
3917 {
3918  Assert(false,
3919  ExcMessage("You are attempting an illegal conversion between "
3920  "iterator/accessor types. The constructor you call "
3921  "only exists to make certain template constructs "
3922  "easier to write as dimension independent code but "
3923  "the conversion is not valid in the current context."));
3924 }
3925 
3926 
3927 
3928 template <int dim, int spacedim>
3929 template <int structdim2, int dim2, int spacedim2>
3932 {
3933  Assert(false,
3934  ExcMessage("You are attempting an illegal conversion between "
3935  "iterator/accessor types. The constructor you call "
3936  "only exists to make certain template constructs "
3937  "easier to write as dimension independent code but "
3938  "the conversion is not valid in the current context."));
3939 }
3940 
3941 
3942 
3943 template <int structdim, int dim, int spacedim>
3944 template <int structdim2, int dim2, int spacedim2>
3947 {
3948  Assert(false,
3949  ExcMessage("You are attempting an illegal conversion between "
3950  "iterator/accessor types. The constructor you call "
3951  "only exists to make certain template constructs "
3952  "easier to write as dimension independent code but "
3953  "the conversion is not valid in the current context."));
3954 }
3955 
3956 
3957 
3958 template <int dim, int spacedim>
3959 template <int structdim2, int dim2, int spacedim2>
3962 {
3963  Assert(false,
3964  ExcMessage("You are attempting an illegal conversion between "
3965  "iterator/accessor types. The constructor you call "
3966  "only exists to make certain template constructs "
3967  "easier to write as dimension independent code but "
3968  "the conversion is not valid in the current context."));
3969 }
3970 
3971 
3972 #ifndef DOXYGEN
3973 
3974 template <>
3975 bool
3977 template <>
3978 bool
3980 template <>
3981 bool
3983 template <>
3984 bool
3986 template <>
3987 bool
3989 template <>
3990 bool
3992 // -------------------------------------------------------------------
3993 
3994 template <>
3995 void
3997 
3998 #endif // DOXYGEN
3999 
4001 
4002 // include more templates in debug and optimized mode
4003 #include "tria_accessor.templates.h"
4004 
4005 #endif
void recursively_set_subdomain_id(const types::subdomain_id new_subdomain_id) const
CellAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
TriaIterator< CellAccessor< dim, spacedim > > parent() const
unsigned int neighbor_face_no(const unsigned int neighbor) const
unsigned int face_iterator_to_index(const TriaIterator< TriaAccessor< dim - 1, dim, spacedim >> &face) const
void set_active_cell_index(const unsigned int active_cell_index) const
unsigned int neighbor_of_neighbor_internal(const unsigned int neighbor) const
TriaIterator< CellAccessor< dim, spacedim > > periodic_neighbor(const unsigned int i) const
types::global_cell_index global_active_cell_index() const
void set_direction_flag(const bool new_direction_flag) const
void recursively_set_material_id(const types::material_id new_material_id) const
void set_level_subdomain_id(const types::subdomain_id new_level_subdomain_id) const
types::subdomain_id level_subdomain_id() const
boost::container::small_vector< TriaIterator< CellAccessor< dim, spacedim > >, GeometryInfo< dim >::max_children_per_cell > child_iterators() const
bool flag_for_face_refinement(const unsigned int face_no, const RefinementCase< dim - 1 > &face_refinement_case=RefinementCase< dim - 1 >::isotropic_refinement) const
unsigned int face_index(const unsigned int i) const
double diameter(const Mapping< dim, spacedim > &mapping) const
::internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
bool is_active() const
bool is_ghost() const
TriaIterator< CellAccessor< dim, spacedim > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
void set_subdomain_id(const types::subdomain_id new_subdomain_id) const
CellAccessor< dim, spacedim > & operator=(CellAccessor< dim, spacedim > &&)=default
bool neighbor_is_coarser(const unsigned int face_no) const
void set_global_level_cell_index(const types::global_cell_index index) const
bool has_periodic_neighbor(const unsigned int i) const
bool active() const
int periodic_neighbor_level(const unsigned int i) const
std::pair< unsigned int, unsigned int > neighbor_of_coarser_neighbor(const unsigned int neighbor) const
~CellAccessor()=default
CellAccessor(const CellAccessor< dim, spacedim > &)=default
void set_coarsen_flag() const
TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > face(const unsigned int i) const
unsigned int neighbor_of_neighbor(const unsigned int face_no) const
void set_material_id(const types::material_id new_material_id) const
bool is_locally_owned() const
TriaIterator< CellAccessor< dim, spacedim > > neighbor(const unsigned int face_no) const
void set_refine_flag(const RefinementCase< dim > ref_case=RefinementCase< dim >::isotropic_refinement) const
RefinementCase< dim > refine_flag_set() const
CellAccessor(CellAccessor< dim, spacedim > &&)=default
bool point_inside_codim(const Point< spacedim_ > &p) const
typename TriaAccessor< dim, dim, spacedim >::AccessorData AccessorData
bool is_locally_owned_on_level() const
bool has_boundary_lines() const
TriaIterator< CellAccessor< dim, spacedim > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
int neighbor_level(const unsigned int face_no) const
int periodic_neighbor_index(const unsigned int i) const
bool periodic_neighbor_is_coarser(const unsigned int i) const
void set_global_active_cell_index(const types::global_cell_index index) const
void clear_coarsen_flag() const
TriaIterator< CellAccessor< dim, spacedim > > child(const unsigned int i) const
void set_neighbor(const unsigned int i, const TriaIterator< CellAccessor< dim, spacedim >> &pointer) const
void set_parent(const unsigned int parent_index)
std::pair< unsigned int, unsigned int > periodic_neighbor_of_coarser_periodic_neighbor(const unsigned face_no) const
CellAccessor(const TriaAccessor< dim, dim, spacedim > &cell_accessor)
bool at_boundary() const
unsigned int active_cell_index() const
static bool is_level_cell()
void clear_refine_flag() const
int neighbor_index(const unsigned int face_no) const
bool point_inside(const Point< spacedim > &p) const
types::subdomain_id subdomain_id() const
bool direction_flag() const
types::material_id material_id() const
bool coarsen_flag_set() const
types::global_cell_index global_level_cell_index() const
bool flag_for_line_refinement(const unsigned int line_no) const
boost::container::small_vector< TriaIterator< TriaAccessor< dim - 1, dim, spacedim > >, GeometryInfo< dim >::faces_per_cell > face_iterators() const
CellId id() const
bool is_artificial() const
TriaIterator< CellAccessor< dim, spacedim > > neighbor_or_periodic_neighbor(const unsigned int i) const
int parent_index() const
unsigned int periodic_neighbor_of_periodic_neighbor(const unsigned int i) const
unsigned int periodic_neighbor_face_no(const unsigned int i) const
CellAccessor< dim, spacedim > & operator=(const CellAccessor< dim, spacedim > &)=delete
Definition: cell_id.h:71
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator quad(const unsigned int i) const
void set_manifold_id(const types::manifold_id) const
InvalidAccessor(const InvalidAccessor &)
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(const unsigned int i) const
InvalidAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
bool operator==(const InvalidAccessor &) const
void operator++() const
bool operator!=(const InvalidAccessor &) const
bool has_children() const
unsigned int user_index() const
Point< spacedim > & vertex(const unsigned int i) const
void copy_from(const InvalidAccessor &)
void set_user_index(const unsigned int p) const
types::manifold_id manifold_id() const
void operator--() const
bool used() const
InvalidAccessor(const OtherAccessor &)
typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
Abstract base class for mapping classes.
Definition: mapping.h:304
bool operator!=(const TriaAccessorBase &) const
static const unsigned int dimension
TriaAccessorBase(const TriaAccessorBase &)
void operator=(const TriaAccessorBase *)=delete
static const unsigned int structure_dimension
typename ::internal::TriaAccessorImplementation::PresentLevelType< structdim, dim >::type present_level
TriaAccessorBase(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *=nullptr)
void copy_from(const TriaAccessorBase &)
IteratorState::IteratorStates state() const
int index() const
bool operator<(const TriaAccessorBase &other) const
TriaAccessorBase & operator=(const TriaAccessorBase &)
static const unsigned int space_dimension
::internal::TriangulationImplementation::TriaObjects & objects() const
int level() const
const Triangulation< dim, spacedim > * tria
const Triangulation< dim, spacedim > & get_triangulation() const
bool operator==(const TriaAccessorBase &) const
static unsigned int line_index(const unsigned int i)
static bool face_orientation(const unsigned int face)
Always return false.
static unsigned int number_of_children()
static RefinementCase< 0 > refinement_case()
static unsigned int max_refinement_depth()
Point< spacedim > center() const
static unsigned int quad_index(const unsigned int i)
TriaAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
static bool face_flip(const unsigned int face)
Always return false.
TriaAccessor(const Triangulation< 1, spacedim > *tria, const VertexKind vertex_kind, const unsigned int vertex_index)
static TriaIterator< TriaAccessor< 0, 1, spacedim > > isotropic_child(const unsigned int)
Return an invalid object.
std_cxx20::ranges::iota_view< unsigned int, unsigned int > line_indices() const
void set_all_boundary_ids(const types::boundary_id)
unsigned int n_lines() const
bool operator!=(const TriaAccessor &) const
TriaAccessor(const Triangulation< 1, spacedim > *tria=nullptr, const int=0, const int=0, const AccessorData *=nullptr)
static unsigned int n_active_descendants()
static int isotropic_child_index(const unsigned int i)
Returns -1.
types::boundary_id boundary_id() const
static typename ::internal::TriangulationImplementation::Iterators< 1, spacedim >::line_iterator line(const unsigned int)
Point< spacedim > & vertex(const unsigned int i=0) const
static bool line_orientation(const unsigned int line)
Always return false.
static typename ::internal::TriangulationImplementation::Iterators< 1, spacedim >::quad_iterator quad(const unsigned int i)
ReferenceCell reference_cell() const
static unsigned int child_iterator_to_index(const TriaIterator< TriaAccessor< 0, 1, spacedim >> &)
Return an invalid unsigned integer.
unsigned int vertex_index(const unsigned int i=0) const
const Triangulation< 1, spacedim > & get_triangulation() const
void copy_from(const TriaAccessor &)
static TriaIterator< TriaAccessor< 0, 1, spacedim > > child(const unsigned int)
Return an invalid object.
static IteratorState::IteratorStates state()
static bool face_rotation(const unsigned int face)
Always return false.
void set_boundary_id(const types::boundary_id)
static int child_index(const unsigned int i)
Returns -1.
const Manifold< 1, spacedim > & get_manifold() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
types::manifold_id manifold_id() const
static unsigned int n_children()
TriaAccessor(const TriaAccessor< structdim2, dim2, spacedim2 > &)
unsigned int n_vertices() const
void set_manifold_id(const types::manifold_id)
bool operator==(const TriaAccessor &) const
const Triangulation< 1, spacedim > * tria
double extent_in_direction(const unsigned int axis) const
TriaAccessor(const Triangulation< dim, spacedim > *tria, const unsigned int vertex_index)
static RefinementCase< 0 > refinement_case()
bool operator!=(const TriaAccessor &) const
static typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(const unsigned int)
void copy_from(const TriaAccessor &)
static typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator quad(const unsigned int i)
TriaAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
static int child_index(const unsigned int i)
Returns -1.
bool operator==(const TriaAccessor &) const
static unsigned int line_index(const unsigned int i)
unsigned int vertex_index(const unsigned int i=0) const
const Triangulation< dim, spacedim > * tria
static int isotropic_child_index(const unsigned int i)
Returns -1.
Point< spacedim > center(const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
TriaAccessor(const Triangulation< dim, spacedim > *tria=nullptr, const int level=0, const int index=0, const AccessorData *=nullptr)
static unsigned int n_children()
static bool face_flip(const unsigned int face)
Always return false.
static TriaIterator< TriaAccessor< 0, dim, spacedim > > isotropic_child(const unsigned int)
Return an invalid object.
static bool line_orientation(const unsigned int line)
Always return false.
static unsigned int child_iterator_to_index(const TriaIterator< TriaAccessor< 0, dim, spacedim >> &)
Return an invalid unsigned integer.
static unsigned int max_refinement_depth()
static unsigned int quad_index(const unsigned int i)
static bool face_rotation(const unsigned int face)
Always return false.
const Triangulation< dim, spacedim > & get_triangulation() const
static bool face_orientation(const unsigned int face)
Always return false.
static unsigned int n_active_descendants()
static TriaIterator< TriaAccessor< 0, dim, spacedim > > child(const unsigned int)
Return an invalid object.
TriaAccessor(const TriaAccessor< structdim2, dim2, spacedim2 > &)
IteratorState::IteratorStates state() const
Point< spacedim > & vertex(const unsigned int i=0) const
static unsigned int number_of_children()
void clear_children() const
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator quad(const unsigned int i) const
bool line_orientation(const unsigned int line) const
void set_boundary_id_internal(const types::boundary_id id) const
void set_user_index(const unsigned int p) const
void set_face_flip(const unsigned int face, const bool flip) const
TriaAccessor(TriaAccessor &&)=default
void clear_user_pointer() const
unsigned int child_iterator_to_index(const TriaIterator< TriaAccessor< structdim, dim, spacedim >> &child) const
unsigned int n_active_descendants() const
void set_face_orientation(const unsigned int face, const bool orientation) const
void set_face_rotation(const unsigned int face, const bool rotation) const
bool face_rotation(const unsigned int face) const
void recursively_set_user_index(const unsigned int p) const
void clear_user_data() const
TriaAccessor(const TriaAccessor< structdim2, dim2, spacedim2 > &)
Point< spacedim > & vertex(const unsigned int i) const
Point< structdim > real_to_unit_cell_affine_approximation(const Point< spacedim > &point) const
void recursively_clear_user_index() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > line_indices() const
unsigned int line_index(const unsigned int i) const
bool face_orientation(const unsigned int face) const
void recursively_set_user_pointer(void *p) const
unsigned int n_lines() const
double extent_in_direction(const unsigned int axis) const
Point< spacedim > intermediate_point(const Point< structdim > &coordinates) const
unsigned int n_vertices() const
bool has_children() const
void recursively_clear_user_flag() const
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(const unsigned int i) const
typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
Point< spacedim > barycenter() const
BoundingBox< spacedim > bounding_box() const
void clear_user_flag() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices() const
TriaAccessor & operator=(TriaAccessor &&)=default
TriaAccessor & operator=(const TriaAccessor &)=delete
unsigned int n_children() const
void recursively_set_user_flag() const
void set_boundary_id(const types::boundary_id) const
bool user_flag_set() const
void set_used_flag() const
types::manifold_id manifold_id() const
void set_user_flag() const
TriaAccessor(const TriaAccessor &)=default
TriaIterator< TriaAccessor< 0, dim, spacedim > > vertex_iterator(const unsigned int i) const
void clear_used_flag() const
std::pair< Point< spacedim >, double > enclosing_ball() const
unsigned int vertex_index(const unsigned int i) const
int isotropic_child_index(const unsigned int i) const
Point< spacedim > center(const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
~TriaAccessor()=default
void set_refinement_case(const RefinementCase< structdim > &ref_case) const
void clear_user_index() const
double minimum_vertex_distance() const
double measure() const
void set_all_boundary_ids(const types::boundary_id) const
TriaIterator< TriaAccessor< structdim, dim, spacedim > > isotropic_child(const unsigned int i) const
void set_bounding_object_indices(const std::initializer_list< int > &new_indices) const
unsigned int number_of_children() const
unsigned int max_refinement_depth() const
void set_line_orientation(const unsigned int line, const bool orientation) const
bool is_translation_of(const TriaIterator< TriaAccessor< structdim, dim, spacedim >> &o) const
unsigned int quad_index(const unsigned int i) const
unsigned int user_index() const
int child_index(const unsigned int i) const
void set_user_pointer(void *p) const
void recursively_clear_user_pointer() const
ReferenceCell reference_cell() const
void * user_pointer() const
TriaIterator< TriaAccessor< structdim, dim, spacedim > > child(const unsigned int i) const
bool face_flip(const unsigned int face) const
TriaAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
const Manifold< dim, spacedim > & get_manifold() const
RefinementCase< structdim > refinement_case() const
TriaAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
bool used() const
void set_children(const unsigned int i, const int index) const
void clear_refinement_case() const
double diameter() const
bool at_boundary() const
types::boundary_id boundary_id() const
unsigned int n_faces() const
#define DEAL_II_DEPRECATED
Definition: config.h:159
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:410
#define DEAL_II_DEPRECATED_EARLY
Definition: config.h:167
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:448
unsigned int level
Definition: grid_out.cc:4590
#define DeclException0(Exception0)
Definition: exceptions.h:470
static ::ExceptionBase & ExcRefineCellNotActive()
static ::ExceptionBase & ExcCantSetChildren(int arg1)
static ::ExceptionBase & ExcCellFlaggedForCoarsening()
static ::ExceptionBase & ExcCellNotActive()
static ::ExceptionBase & ExcCellHasNoParent()
static ::ExceptionBase & ExcNeighborIsNotCoarser()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcCellNotUsed()
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcCellFlaggedForRefinement()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcDereferenceInvalidObject(AccessorType arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
static ::ExceptionBase & ExcSetOnlyEvenChildren(int arg1)
static ::ExceptionBase & ExcCellHasNoChildren()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcFacesHaveNoLevel()
static ::ExceptionBase & ExcNoPeriodicNeighbor()
static ::ExceptionBase & ExcNeighborIsCoarser()
void set_all_manifold_ids(const types::manifold_id) const
void set_all_manifold_ids(const types::manifold_id)
void set_manifold_id(const types::manifold_id) const
@ past_the_end
Iterator reached end of container.
@ valid
Iterator points to a valid object.
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
boost::integer_range< IncrementableType > iota_view
Definition: iota_view.h:46
unsigned int manifold_id
Definition: types.h:141
unsigned int subdomain_id
Definition: types.h:43
unsigned int global_cell_index
Definition: types.h:105
unsigned int material_id
Definition: types.h:152
unsigned int boundary_id
Definition: types.h:129