Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_tria_accessor_h
17 #define dealii_tria_accessor_h
18 
19 
20 #include <deal.II/base/bounding_box.h>
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/geometry_info.h>
24 #include <deal.II/base/point.h>
25 #include <deal.II/grid/tria_iterator_base.h>
26 #include <deal.II/grid/tria_iterator_selector.h>
27 #include <deal.II/grid/cell_id.h>
28 
29 #include <utility>
30 
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 template <int dim, int spacedim> class Triangulation;
35 template <typename Accessor> class TriaRawIterator;
36 template <typename Accessor> class TriaIterator;
37 template <typename Accessor> class TriaActiveIterator;
38 
39 template <int dim, int spacedim> class Boundary;
40 template <int dim, int spacedim> class Manifold;
41 
42 
43 namespace internal
44 {
45  namespace TriangulationImplementation
46  {
47  template <int dim> class TriaObject;
48  template <typename G> class TriaObjects;
49  struct Implementation;
50  }
51 
52  namespace TriaAccessorImplementation
53  {
54  struct Implementation;
55 
61  template <int structdim, int dim> struct PresentLevelType
62  {
63  struct type
64  {
68  type () = default;
69 
73  type (const int level)
74  {
75  Assert (level == 0, ExcInternalError());
76  (void)level; // removes -Wunused-parameter warning in optimized mode
77  }
78 
82  operator int () const
83  {
84  return 0;
85  }
86 
87  void operator ++ () const
88  {
89  Assert (false, ExcInternalError());
90  }
91 
92  void operator -- () const
93  {
94  Assert (false, ExcInternalError());
95  }
96  };
97  };
98 
99 
105  template <int dim> struct PresentLevelType<dim,dim>
106  {
107  typedef int type;
108  };
109 
110  }
111 }
112 template <int structdim, int dim, int spacedim> class TriaAccessor;
113 template <int dim, int spacedim> class TriaAccessor<0, dim, spacedim>;
114 template <int spacedim> class TriaAccessor<0, 1, spacedim>;
115 
116 // note: the file tria_accessor.templates.h is included at the end of
117 // this file. this includes a lot of templates. originally, this was
118 // only done in debug mode, but led to cyclic reduction problems and
119 // so is now on by default.
120 
121 
126 {
131  "The operation you are attempting can only be performed for "
132  "(cell, face, or edge) iterators that point to valid "
133  "objects. These objects need not necessarily be active, "
134  "i.e., have no children, but they need to be part of a "
135  "triangulation. (The objects pointed to by an iterator "
136  "may -- after coarsening -- also be objects that used "
137  "to be part of a triangulation, but are now no longer "
138  "used. Their memory location may have been retained "
139  "for re-use upon the next mesh refinement, but is "
140  "currently unused.)");
151  "The operation you are attempting can only be performed for "
152  "(cell, face, or edge) iterators that point to 'active' "
153  "objects. 'Active' objects are those that do not have "
154  "children (in the case of cells), or that are part of "
155  "an active cell (in the case of faces or edges). However, "
156  "the object on which you are trying the current "
157  "operation is not 'active' in this sense.");
164  "The operation you are attempting can only be performed for "
165  "(cell, face, or edge) iterators that have children, "
166  "but the object on which you are trying the current "
167  "operation does not have any.");
175  "The operation you are attempting can only be performed for "
176  "(cell, face, or edge) iterators that have a parent object, "
177  "but the object on which you are trying the current "
178  "operation does not have one -- i.e., it is on the "
179  "coarsest level of the triangulation.");
184  int,
185  << "You can only set the child index if the cell does not "
186  << "currently have children registered; or you can clear it. "
187  << "The given index was " << arg1 << " (-1 means: clear children).");
191  template <typename AccessorType>
193  AccessorType,
194  << "You tried to dereference an iterator for which this "
195  << "is not possible. More information on this iterator: "
196  << "index=" << arg1.index()
197  << ", state="
198  << (arg1.state() == IteratorState::valid ? "valid" :
199  (arg1.state() == IteratorState::past_the_end ?
200  "past_the_end" : "invalid")));
205  "Iterators can only be compared if they point to the same "
206  "triangulation, or if neither of them are associated "
207  "with a triangulation.");
208 //TODO: Write documentation!
213 //TODO: Write documentation!
233 //TODO: Write documentation!
238  int,
239  << "You can only set the child index of an even numbered child."
240  << "The number of the child given was " << arg1 << ".");
241 }
242 
243 
270 template <int structdim, int dim, int spacedim=dim>
271 class TriaAccessorBase
272 {
273 public:
279  static const unsigned int space_dimension = spacedim;
280 
286  static const unsigned int dimension = dim;
287 
293  static const unsigned int structure_dimension = structdim;
294 
295 protected:
301  typedef void AccessorData;
302 
306  TriaAccessorBase (const Triangulation<dim,spacedim> *parent = nullptr,
307  const int level = -1,
308  const int index = -1,
309  const AccessorData * = nullptr);
310 
315 
324  void operator = (const TriaAccessorBase *) = delete;
325 
333  void copy_from (const TriaAccessorBase &);
334 
339 
346  bool operator < (const TriaAccessorBase &other) const;
347 
348 protected:
349 
353  bool operator == (const TriaAccessorBase &) const;
354 
358  bool operator != (const TriaAccessorBase &) const;
359 
373  void operator ++ ();
374 
382  void operator -- ();
391  objects () const;
392 
393 public:
399  typedef void *LocalData;
400 
424  int level () const;
425 
452  int index () const;
453 
459 
465 
469 protected:
474  typename ::internal::TriaAccessorImplementation::PresentLevelType<structdim,dim>::type present_level;
475 
481 
486 
487 private:
488 
489  template <typename Accessor> friend class TriaRawIterator;
490  template <typename Accessor> friend class TriaIterator;
491  template <typename Accessor> friend class TriaActiveIterator;
492 };
493 
494 
495 
517 template <int structdim, int dim, int spacedim=dim>
518 class InvalidAccessor : public TriaAccessorBase<structdim,dim,spacedim>
519 {
520 public:
525 
533  InvalidAccessor (const Triangulation<dim,spacedim> *parent = nullptr,
534  const int level = -1,
535  const int index = -1,
536  const AccessorData *local_data = nullptr);
537 
546 
551  template <typename OtherAccessor>
552  InvalidAccessor (const OtherAccessor &);
553 
557  void copy_from (const InvalidAccessor &);
558 
562  bool operator == (const InvalidAccessor &) const;
563  bool operator != (const InvalidAccessor &) const;
564 
568  void operator ++ () const;
569  void operator -- () const;
570 
575  bool used () const;
576 
581  bool has_children () const;
582 
587 
591  Point<spacedim> &vertex (const unsigned int i) const;
592 
596  typename ::internal::TriangulationImplementation::Iterators<dim,spacedim>::line_iterator
597  line (const unsigned int i) const;
598 
602  typename ::internal::TriangulationImplementation::Iterators<dim,spacedim>::quad_iterator
603  quad (const unsigned int i) const;
604 };
605 
606 
607 
626 template <int structdim, int dim, int spacedim>
627 class TriaAccessor : public TriaAccessorBase<structdim, dim, spacedim>
628 {
629 public:
633  typedef
636 
640  TriaAccessor (const Triangulation<dim,spacedim> *parent = nullptr,
641  const int level = -1,
642  const int index = -1,
643  const AccessorData *local_data = nullptr);
644 
657  template <int structdim2, int dim2, int spacedim2>
659 
664  template <int structdim2, int dim2, int spacedim2>
666 
675  void operator = (const TriaAccessor &) = delete;
676 
683  bool used () const;
684 
697  vertex_iterator (const unsigned int i) const;
698 
714  unsigned int vertex_index (const unsigned int i) const;
715 
753  Point<spacedim> &vertex (const unsigned int i) const;
754 
758  typename ::internal::TriangulationImplementation::Iterators<dim,spacedim>::line_iterator
759  line (const unsigned int i) const;
760 
767  unsigned int line_index (const unsigned int i) const;
768 
772  typename ::internal::TriangulationImplementation::Iterators<dim,spacedim>::quad_iterator
773  quad (const unsigned int i) const;
774 
781  unsigned int quad_index (const unsigned int i) const;
804  bool face_orientation (const unsigned int face) const;
805 
815  bool face_flip (const unsigned int face) const;
816 
826  bool face_rotation (const unsigned int face) const;
827 
838  bool line_orientation (const unsigned int line) const;
853  bool has_children () const;
854 
859  unsigned int n_children() const;
860 
874  unsigned int number_of_children () const;
875 
889  unsigned int max_refinement_depth () const;
890 
895  child (const unsigned int i) const;
896 
906  isotropic_child (const unsigned int i) const;
907 
912 
918  int child_index (const unsigned int i) const;
919 
925  int isotropic_child_index (const unsigned int i) const;
948 
978  void set_boundary_id (const types::boundary_id) const;
979 
1010  void set_all_boundary_ids (const types::boundary_id) const;
1011 
1019  bool at_boundary () const;
1020 
1029  DEAL_II_DEPRECATED const Boundary<dim,spacedim> &get_boundary () const;
1030 
1040  const Manifold<dim,spacedim> &get_manifold () const;
1041 
1063 
1081  void set_manifold_id (const types::manifold_id) const;
1082 
1096  void set_all_manifold_ids (const types::manifold_id) const;
1097 
1114  bool user_flag_set () const;
1115 
1121  void set_user_flag () const;
1122 
1128  void clear_user_flag () const;
1129 
1135  void recursively_set_user_flag () const;
1136 
1142  void recursively_clear_user_flag () const;
1143 
1149  void clear_user_data () const;
1150 
1162  void set_user_pointer (void *p) const;
1163 
1169  void clear_user_pointer () const;
1170 
1186  void *user_pointer () const;
1187 
1209  void recursively_set_user_pointer (void *p) const;
1210 
1217  void recursively_clear_user_pointer () const;
1218 
1228  void set_user_index (const unsigned int p) const;
1229 
1235  void clear_user_index () const;
1236 
1248  unsigned int user_index () const;
1249 
1267  void recursively_set_user_index (const unsigned int p) const;
1268 
1277  void recursively_clear_user_index () const;
1296  double diameter () const;
1297 
1324  std::pair<Point<spacedim>,double> enclosing_ball () const;
1325 
1330 
1340  double extent_in_direction (const unsigned int axis) const;
1341 
1345  double minimum_vertex_distance () const;
1346 
1360  Point<spacedim> intermediate_point(const Point<structdim> &coordinates) const;
1361 
1386 
1421  Point<spacedim> center (const bool respect_manifold=false,
1422  const bool interpolate_from_surrounding=false) const;
1423 
1441  Point<spacedim> barycenter () const;
1442 
1468  double measure () const;
1469 
1484  bool
1486 
1492 private:
1496  void set_boundary_id_internal(const types::boundary_id id) const;
1497 
1502  void set (const ::internal::TriangulationImplementation::TriaObject<structdim> &o) const;
1503 
1511  void set_line_orientation (const unsigned int line,
1512  const bool orientation) const;
1513 
1524  void set_face_orientation (const unsigned int face,
1525  const bool orientation) const;
1526 
1533  void set_face_flip (const unsigned int face,
1534  const bool flip) const;
1535 
1542  void set_face_rotation (const unsigned int face,
1543  const bool rotation) const;
1544 
1548  void set_used_flag () const;
1549 
1553  void clear_used_flag () const;
1554 
1563  void set_refinement_case (const RefinementCase<structdim> &ref_case) const;
1564 
1572  void clear_refinement_case () const;
1573 
1580  void set_children (const unsigned int i, const int index) const;
1581 
1586  void clear_children () const;
1587 
1588 private:
1589 
1590  template <int, int> friend class Triangulation;
1591 
1592  friend struct ::internal::TriangulationImplementation::Implementation;
1593  friend struct ::internal::TriaAccessorImplementation::Implementation;
1594 };
1595 
1596 
1597 
1598 
1599 
1600 
1619 template <int dim, int spacedim>
1620 class TriaAccessor<0, dim, spacedim>
1621 {
1622 public:
1628  static const unsigned int space_dimension = spacedim;
1629 
1635  static const unsigned int dimension = dim;
1636 
1642  static const unsigned int structure_dimension = 0;
1643 
1647  typedef void AccessorData;
1648 
1654  const unsigned int vertex_index);
1655 
1661  TriaAccessor (const Triangulation<dim,spacedim> *tria = nullptr,
1662  const int level = 0,
1663  const int index = 0,
1664  const AccessorData * = nullptr);
1665 
1669  template <int structdim2, int dim2, int spacedim2>
1671 
1675  template <int structdim2, int dim2, int spacedim2>
1677 
1682 
1687  static int level ();
1688 
1693  int index () const;
1694 
1704  void operator ++ ();
1705 
1709  void operator -- ();
1713  bool operator == (const TriaAccessor &) const;
1714 
1718  bool operator != (const TriaAccessor &) const;
1719 
1747  unsigned int vertex_index (const unsigned int i = 0) const;
1748 
1754  Point<spacedim> &vertex (const unsigned int i = 0) const;
1755 
1760  typename ::internal::TriangulationImplementation::Iterators<dim,spacedim>::line_iterator
1761  static line (const unsigned int);
1762 
1766  static unsigned int line_index (const unsigned int i);
1767 
1771  static
1772  typename ::internal::TriangulationImplementation::Iterators<dim,spacedim>::quad_iterator
1773  quad (const unsigned int i);
1774 
1778  static unsigned int quad_index (const unsigned int i);
1779 
1795  double diameter () const;
1796 
1804  double extent_in_direction (const unsigned int axis) const;
1805 
1813  Point<spacedim> center (const bool respect_manifold=false,
1814  const bool interpolate_from_surrounding=false) const;
1815 
1824  double measure () const;
1839  static bool face_orientation (const unsigned int face);
1840 
1844  static bool face_flip (const unsigned int face);
1845 
1849  static bool face_rotation (const unsigned int face);
1850 
1854  static bool line_orientation (const unsigned int line);
1855 
1870  static bool has_children ();
1871 
1876  static unsigned int n_children();
1877 
1882  static unsigned int number_of_children ();
1883 
1887  static unsigned int max_refinement_depth ();
1888 
1892  static
1894  child (const unsigned int);
1895 
1899  static
1901  isotropic_child (const unsigned int);
1902 
1906  static
1908 
1912  static
1913  int child_index (const unsigned int i);
1914 
1918  static
1919  int isotropic_child_index (const unsigned int i);
1927  bool used () const;
1928 
1929 protected:
1937  void copy_from (const TriaAccessor &);
1938 
1943 
1947  unsigned int global_vertex_index;
1948 
1949 private:
1950 
1951  template <typename Accessor> friend class TriaRawIterator;
1952  template <typename Accessor> friend class TriaIterator;
1953  template <typename Accessor> friend class TriaActiveIterator;
1954 };
1955 
1956 
1957 
1974 template <int spacedim>
1975 class TriaAccessor<0, 1, spacedim>
1976 {
1977 public:
1983  static const unsigned int space_dimension = spacedim;
1984 
1990  static const unsigned int dimension = 1;
1991 
1997  static const unsigned int structure_dimension = 0;
1998 
2002  typedef void AccessorData;
2003 
2009  {
2021  right_vertex
2022  };
2023 
2036  const VertexKind vertex_kind,
2037  const unsigned int vertex_index);
2038 
2044  TriaAccessor (const Triangulation<1,spacedim> *tria = nullptr,
2045  const int = 0,
2046  const int = 0,
2047  const AccessorData * = nullptr);
2048 
2052  template <int structdim2, int dim2, int spacedim2>
2054 
2058  template <int structdim2, int dim2, int spacedim2>
2060 
2065  void copy_from (const TriaAccessor &);
2066 
2073 
2078  static int level ();
2079 
2084  int index () const;
2085 
2096  void operator ++ () const;
2097 
2102  void operator -- () const;
2106  bool operator == (const TriaAccessor &) const;
2107 
2111  bool operator != (const TriaAccessor &) const;
2112 
2139  unsigned int vertex_index (const unsigned int i = 0) const;
2140 
2146  Point<spacedim> &vertex (const unsigned int i = 0) const;
2147 
2152  Point<spacedim> center () const;
2153 
2158  typename ::internal::TriangulationImplementation::Iterators<1,spacedim>::line_iterator
2159  static line (const unsigned int);
2160 
2167  static unsigned int line_index (const unsigned int i);
2168 
2172  static
2173  typename ::internal::TriangulationImplementation::Iterators<1,spacedim>::quad_iterator
2174  quad (const unsigned int i);
2175 
2182  static unsigned int quad_index (const unsigned int i);
2183 
2193  bool at_boundary () const;
2194 
2210 
2214  const Manifold<1,spacedim> &get_manifold () const;
2215 
2223 
2224 
2235  static bool face_orientation (const unsigned int face);
2236 
2240  static bool face_flip (const unsigned int face);
2241 
2245  static bool face_rotation (const unsigned int face);
2246 
2250  static bool line_orientation (const unsigned int line);
2251 
2266  static bool has_children ();
2267 
2272  static unsigned int n_children();
2273 
2278  static unsigned int number_of_children ();
2279 
2283  static unsigned int max_refinement_depth ();
2284 
2288  static
2290  child (const unsigned int);
2291 
2295  static
2297  isotropic_child (const unsigned int);
2298 
2302  static
2304 
2308  static
2309  int child_index (const unsigned int i);
2310 
2314  static
2315  int isotropic_child_index (const unsigned int i);
2348  void
2350 
2357  void
2359 
2371  void
2373 
2385  void
2394  bool used () const;
2395 
2396 protected:
2401 
2407 
2411  unsigned int global_vertex_index;
2412 };
2413 
2414 
2415 
2416 
2433 template <int dim, int spacedim=dim>
2434 class CellAccessor : public TriaAccessor<dim,dim,spacedim>
2435 {
2436 public:
2441 
2446 
2458  const int level = -1,
2459  const int index = -1,
2460  const AccessorData *local_data = nullptr);
2461 
2465  CellAccessor (const TriaAccessor<dim,dim,spacedim> &cell_accessor);
2466 
2479  template <int structdim2, int dim2, int spacedim2>
2481 
2486  template <int structdim2, int dim2, int spacedim2>
2488 
2497  void operator = (const CellAccessor<dim, spacedim> &) = delete;
2498 
2515  child (const unsigned int i) const;
2516 
2520  TriaIterator<TriaAccessor<dim-1,dim,spacedim> >
2521  face (const unsigned int i) const;
2522 
2532  unsigned int
2533  face_index (const unsigned int i) const;
2534 
2583  neighbor_child_on_subface (const unsigned int face_no,
2584  const unsigned int subface_no) const;
2585 
2608  neighbor (const unsigned int i) const;
2609 
2614  int neighbor_index (const unsigned int i) const;
2615 
2620  int neighbor_level (const unsigned int i) const;
2621 
2633  unsigned int neighbor_of_neighbor (const unsigned int neighbor) const;
2634 
2645  bool neighbor_is_coarser (const unsigned int neighbor) const;
2646 
2661  std::pair<unsigned int, unsigned int>
2662  neighbor_of_coarser_neighbor (const unsigned int neighbor) const;
2663 
2670  unsigned int neighbor_face_no (const unsigned int neighbor) const;
2671 
2675  static bool is_level_cell();
2676 
2690  bool has_periodic_neighbor(const unsigned int i) const;
2691 
2709  periodic_neighbor (const unsigned int i) const;
2710 
2719  neighbor_or_periodic_neighbor (const unsigned int i) const;
2720 
2736  periodic_neighbor_child_on_subface (const unsigned int face_no,
2737  const unsigned int subface_no) const;
2738 
2748  std::pair<unsigned int, unsigned int>
2749  periodic_neighbor_of_coarser_periodic_neighbor (const unsigned face_no) const;
2750 
2756  int
2757  periodic_neighbor_index (const unsigned int i) const;
2758 
2764  int
2765  periodic_neighbor_level (const unsigned int i) const;
2766 
2781  unsigned int
2782  periodic_neighbor_of_periodic_neighbor (const unsigned int i) const;
2783 
2789  unsigned int
2790  periodic_neighbor_face_no (const unsigned int i) const;
2791 
2798  bool
2799  periodic_neighbor_is_coarser (const unsigned int i) const;
2800 
2817  bool at_boundary (const unsigned int i) const;
2818 
2827  bool at_boundary () const;
2828 
2836  bool has_boundary_lines () const;
2863 
2882 
2886  void clear_refine_flag () const;
2887 
2895  bool flag_for_face_refinement (const unsigned int face_no,
2896  const RefinementCase<dim-1> &face_refinement_case=RefinementCase<dim-1>::isotropic_refinement) const;
2897 
2903  bool flag_for_line_refinement (const unsigned int line_no) const;
2904 
2913  ::internal::SubfaceCase<dim> subface_case(const unsigned int face_no) const;
2914 
2918  bool coarsen_flag_set () const;
2919 
2924  void set_coarsen_flag () const;
2925 
2929  void clear_coarsen_flag () const;
2953 
2965  void set_material_id (const types::material_id new_material_id) const;
2966 
2975  void recursively_set_material_id (const types::material_id new_material_id) const;
3002 
3018  void set_subdomain_id (const types::subdomain_id new_subdomain_id) const;
3019 
3025 
3030  void set_level_subdomain_id (const types::subdomain_id new_level_subdomain_id) const;
3031 
3032 
3048  void recursively_set_subdomain_id (const types::subdomain_id new_subdomain_id) const;
3066  bool direction_flag () const;
3067 
3093  unsigned int active_cell_index () const;
3094 
3102  int parent_index () const;
3103 
3110  parent () const;
3111 
3131  bool active () const;
3132 
3152  bool is_locally_owned () const;
3153 
3158  bool is_locally_owned_on_level () const;
3159 
3183  bool is_ghost () const;
3184 
3211  bool is_artificial () const;
3212 
3226  bool point_inside (const Point<spacedim> &p) const;
3227 
3236  void set_neighbor (const unsigned int i,
3237  const TriaIterator<CellAccessor<dim, spacedim> > &pointer) const;
3238 
3252  CellId id() const;
3253 
3271 
3272 protected:
3288  unsigned int neighbor_of_neighbor_internal (const unsigned int neighbor) const;
3289 
3295  template <int dim_,int spacedim_ >
3296  bool point_inside_codim(const Point<spacedim_> &p) const;
3297 
3298 
3299 
3300 private:
3305  void set_active_cell_index (const unsigned int active_cell_index);
3306 
3310  void set_parent (const unsigned int parent_index);
3311 
3318  void set_direction_flag (const bool new_direction_flag) const;
3319 
3320  template <int, int> friend class Triangulation;
3321 
3322  friend struct ::internal::TriangulationImplementation::Implementation;
3323 };
3324 
3325 
3326 
3327 /* -------------- declaration of explicit
3328  specializations and general templates ------------- */
3329 
3330 
3331 template <int structdim, int dim, int spacedim>
3332 template <typename OtherAccessor>
3334 InvalidAccessor (const OtherAccessor &)
3335 {
3336  Assert (false,
3337  ExcMessage ("You are attempting an illegal conversion between "
3338  "iterator/accessor types. The constructor you call "
3339  "only exists to make certain template constructs "
3340  "easier to write as dimension independent code but "
3341  "the conversion is not valid in the current context."));
3342 }
3343 
3344 
3345 
3346 template <int structdim, int dim, int spacedim>
3347 template <int structdim2, int dim2, int spacedim2>
3350 {
3351  Assert (false,
3352  ExcMessage ("You are attempting an illegal conversion between "
3353  "iterator/accessor types. The constructor you call "
3354  "only exists to make certain template constructs "
3355  "easier to write as dimension independent code but "
3356  "the conversion is not valid in the current context."));
3357 }
3358 
3359 
3360 
3361 template <int dim, int spacedim>
3362 template <int structdim2, int dim2, int spacedim2>
3365 {
3366  Assert (false,
3367  ExcMessage ("You are attempting an illegal conversion between "
3368  "iterator/accessor types. The constructor you call "
3369  "only exists to make certain template constructs "
3370  "easier to write as dimension independent code but "
3371  "the conversion is not valid in the current context."));
3372 }
3373 
3374 
3375 
3376 template <int structdim, int dim, int spacedim>
3377 template <int structdim2, int dim2, int spacedim2>
3380 {
3381  Assert (false,
3382  ExcMessage ("You are attempting an illegal conversion between "
3383  "iterator/accessor types. The constructor you call "
3384  "only exists to make certain template constructs "
3385  "easier to write as dimension independent code but "
3386  "the conversion is not valid in the current context."));
3387 }
3388 
3389 
3390 
3391 template <int dim, int spacedim>
3392 template <int structdim2, int dim2, int spacedim2>
3395 {
3396  Assert (false,
3397  ExcMessage ("You are attempting an illegal conversion between "
3398  "iterator/accessor types. The constructor you call "
3399  "only exists to make certain template constructs "
3400  "easier to write as dimension independent code but "
3401  "the conversion is not valid in the current context."));
3402 }
3403 
3404 
3405 #ifndef DOXYGEN
3406 
3407 template <> bool CellAccessor<1,1>::point_inside (const Point<1> &) const;
3408 template <> bool CellAccessor<2,2>::point_inside (const Point<2> &) const;
3409 template <> bool CellAccessor<3,3>::point_inside (const Point<3> &) const;
3410 template <> bool CellAccessor<1,2>::point_inside (const Point<2> &) const;
3411 template <> bool CellAccessor<1,3>::point_inside (const Point<3> &) const;
3412 template <> bool CellAccessor<2,3>::point_inside (const Point<3> &) const;
3413 // -------------------------------------------------------------------
3414 
3415 template <> void TriaAccessor<3,3,3>::set_all_manifold_ids (const types::manifold_id) const;
3416 
3417 #endif // DOXYGEN
3418 
3419 DEAL_II_NAMESPACE_CLOSE
3420 
3421 // include more templates in debug and optimized mode
3422 # include "tria_accessor.templates.h"
3423 
3424 #endif
unsigned int periodic_neighbor_of_periodic_neighbor(const unsigned int i) const
void clear_used_flag() const
TriaAccessor< dim, dim, spacedim >::AccessorData AccessorData
std::pair< Point< spacedim >, double > enclosing_ball() const
static ::ExceptionBase & ExcCellHasNoChildren()
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator quad(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 point_inside_codim(const Point< spacedim_ > &p) const
::internal::TriangulationImplementation::TriaObjects<::internal::TriangulationImplementation::TriaObject< structdim > > & objects() const
void set_face_rotation(const unsigned int face, const bool rotation) const
TriaIterator< CellAccessor< dim, spacedim > > neighbor_or_periodic_neighbor(const unsigned int i) const
static ::ExceptionBase & ExcNeighborIsNotCoarser()
void set_active_cell_index(const unsigned int active_cell_index)
void recursively_clear_user_pointer() const
unsigned int active_cell_index() const
bool operator==(const InvalidAccessor &) const
static ::ExceptionBase & ExcNeighborIsCoarser()
static ::ExceptionBase & ExcRefineCellNotActive()
unsigned int neighbor_of_neighbor_internal(const unsigned int neighbor) const
void set_neighbor(const unsigned int i, const TriaIterator< CellAccessor< dim, spacedim > > &pointer) const
static ::ExceptionBase & ExcDereferenceInvalidObject(AccessorType arg1)
unsigned int user_index() const
static bool is_level_cell()
int isotropic_child_index(const unsigned int i) const
const Triangulation< dim, spacedim > * tria
void * user_pointer() const
types::material_id material_id() const
bool at_boundary() const
TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
bool flag_for_line_refinement(const unsigned int line_no) const
void set_line_orientation(const unsigned int line, const bool orientation) const
unsigned int boundary_id
Definition: types.h:110
TriaIterator< TriaAccessor< 0, dim, spacedim > > vertex_iterator(const unsigned int i) const
TriaIterator< TriaAccessor< dim-1, dim, spacedim > > face(const unsigned int i) const
void copy_from(const TriaAccessorBase &)
double extent_in_direction(const unsigned int axis) const
bool operator!=(const TriaAccessorBase &) const
unsigned int n_children() const
const Boundary< dim, spacedim > & get_boundary() const
void set_all_manifold_ids(const types::manifold_id) const
const Triangulation< 1, spacedim > * tria
void set_children(const unsigned int i, const int index) const
void clear_user_pointer() const
unsigned int line_index(const unsigned int i) const
unsigned int material_id
Definition: types.h:133
static ::ExceptionBase & ExcFacesHaveNoLevel()
void operator=(const CellAccessor< dim, spacedim > &)=delete
unsigned int face_index(const unsigned int i) const
TriaIterator< CellAccessor< dim, spacedim > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(const unsigned int i) const
bool neighbor_is_coarser(const unsigned int neighbor) const
bool has_children() const
void operator++() const
int level() const
unsigned int vertex_index(const unsigned int i) const
types::manifold_id manifold_id() const
bool is_locally_owned() const
bool is_translation_of(const TriaIterator< TriaAccessor< structdim, dim, spacedim > > &o) const
bool face_rotation(const unsigned int face) const
static ::ExceptionBase & ExcCellNotUsed()
void copy_from(const InvalidAccessor &)
bool at_boundary() const
void clear_user_data() const
TriaIterator< CellAccessor< dim, spacedim > > parent() const
static const unsigned int space_dimension
CellId id() const
void set_user_pointer(void *p) const
void set_boundary_id_internal(const types::boundary_id id) const
bool coarsen_flag_set() const
TriaIterator< CellAccessor< dim, spacedim > > neighbor(const unsigned int i) const
Definition: tria.h:45
TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
static ::ExceptionBase & ExcMessage(std::string arg1)
bool has_periodic_neighbor(const unsigned int i) const
void recursively_set_user_flag() const
bool operator==(const TriaAccessorBase &) const
void set_coarsen_flag() const
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(const unsigned int i) const
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:346
void clear_user_flag() const
types::manifold_id manifold_id() const
bool active() const
bool is_locally_owned_on_level() const
bool periodic_neighbor_is_coarser(const unsigned int i) const
static ::ExceptionBase & ExcCellFlaggedForCoarsening()
static ::ExceptionBase & ExcSetOnlyEvenChildren(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1142
std::pair< unsigned int, unsigned int > periodic_neighbor_of_coarser_periodic_neighbor(const unsigned face_no) const
Point< spacedim > & vertex(const unsigned int i) const
BoundingBox< spacedim > bounding_box() const
bool line_orientation(const unsigned int line) const
void set_material_id(const types::material_id new_material_id) const
bool is_ghost() const
Point< spacedim > & vertex(const unsigned int i) const
const Triangulation< dim, spacedim > * tria
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:335
#define DeclException0(Exception0)
Definition: exceptions.h:323
void set_face_orientation(const unsigned int face, const bool orientation) const
unsigned int quad_index(const unsigned int i) const
bool operator<(const TriaAccessorBase &other) const
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator quad(const unsigned int i) const
RefinementCase< dim > refine_flag_set() const
int index() const
static ::ExceptionBase & ExcCellNotActive()
void operator=(const TriaAccessorBase *)=delete
void operator=(const TriaAccessor &)=delete
IteratorState::IteratorStates state() const
void set_refine_flag(const RefinementCase< dim > ref_case=RefinementCase< dim >::isotropic_refinement) const
static ::ExceptionBase & ExcNoPeriodicNeighbor()
static ::ExceptionBase & ExcCantSetChildren(int arg1)
void clear_refine_flag() const
void recursively_set_material_id(const types::material_id new_material_id) const
typename ::internal::TriaAccessorImplementation::PresentLevelType< structdim, dim >::type present_level
unsigned int subdomain_id
Definition: types.h:42
void clear_children() const
void recursively_clear_user_index() const
TriaIterator< CellAccessor< dim, spacedim > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
void set_user_index(const unsigned int p) const
double minimum_vertex_distance() const
void clear_user_index() const
bool has_children() const
void recursively_clear_user_flag() const
Definition: cell_id.h:61
CellAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
Point< spacedim > center(const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
int parent_index() const
TriaAccessorBase(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *=nullptr)
static ::ExceptionBase & ExcCellHasNoParent()
void set_all_boundary_ids(const types::boundary_id) const
Point< structdim > real_to_unit_cell_affine_approximation(const Point< spacedim > &point) const
bool is_artificial() const
types::boundary_id boundary_id() const
double diameter() const
int periodic_neighbor_level(const unsigned int i) const
void recursively_set_user_pointer(void *p) const
std::pair< unsigned int, unsigned int > neighbor_of_coarser_neighbor(const unsigned int neighbor) const
bool flag_for_face_refinement(const unsigned int face_no, const RefinementCase< dim-1 > &face_refinement_case=RefinementCase< dim-1 >::isotropic_refinement) const
void set_manifold_id(const types::manifold_id) const
unsigned int manifold_id
Definition: types.h:122
void recursively_set_subdomain_id(const types::subdomain_id new_subdomain_id) const
void set_boundary_id(const types::boundary_id) const
bool face_flip(const unsigned int face) const
::internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
int neighbor_level(const unsigned int i) const
bool point_inside(const Point< spacedim > &p) const
Point< spacedim > barycenter() const
bool used() const
void set_direction_flag(const bool new_direction_flag) const
static const unsigned int structure_dimension
unsigned int periodic_neighbor_face_no(const unsigned int i) const
int periodic_neighbor_index(const unsigned int i) const
TriaIterator< CellAccessor< dim, spacedim > > periodic_neighbor(const unsigned int i) const
unsigned int number_of_children() const
Triangulation< dim, spacedim > Container
void clear_refinement_case() const
int child_index(const unsigned int i) const
void set_parent(const unsigned int parent_index)
void clear_coarsen_flag() const
bool used() const
bool face_orientation(const unsigned int face) const
Point< spacedim > intermediate_point(const Point< structdim > &coordinates) const
TriaAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
bool direction_flag() const
int neighbor_index(const unsigned int i) const
const Triangulation< dim, spacedim > & get_triangulation() const
Iterator reached end of container.
Iterator points to a valid object.
void set_level_subdomain_id(const types::subdomain_id new_level_subdomain_id) const
void recursively_set_user_index(const unsigned int p) const
TriaIterator< TriaAccessor< structdim, dim, spacedim > > isotropic_child(const unsigned int i) const
void set_subdomain_id(const types::subdomain_id new_subdomain_id) const
types::subdomain_id level_subdomain_id() const
void set_user_flag() const
bool user_flag_set() const
void set_face_flip(const unsigned int face, const bool flip) const
unsigned int neighbor_face_no(const unsigned int neighbor) const
types::subdomain_id subdomain_id() const
static ::ExceptionBase & ExcCellFlaggedForRefinement()
static const unsigned int dimension
double measure() const
const Manifold< dim, spacedim > & get_manifold() const
void set_refinement_case(const RefinementCase< structdim > &ref_case) const
RefinementCase< structdim > refinement_case() const
unsigned int max_refinement_depth() const
unsigned int neighbor_of_neighbor(const unsigned int neighbor) const
void set_used_flag() const
TriaIterator< CellAccessor< dim, spacedim > > child(const unsigned int i) const
bool has_boundary_lines() const
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcCantCompareIterators()
TriaIterator< TriaAccessor< structdim, dim, spacedim > > child(const unsigned int i) const