Reference documentation for deal.II version 8.1.0
tria_accessor.h
1 // ---------------------------------------------------------------------
2 // @f$Id: tria_accessor.h 30264 2013-08-09 12:37:28Z maier @f$
3 //
4 // Copyright (C) 1998 - 2013 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 #ifndef __deal2__tria_accessor_h
18 #define __deal2__tria_accessor_h
19 
20 
21 #include <deal.II/base/config.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 
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 
41 
42 namespace internal
43 {
44  namespace Triangulation
45  {
46  template <int dim> class TriaObject;
47  template <typename G> class TriaObjects;
48  struct Implementation;
49  }
50 
51  namespace TriaAccessor
52  {
53  struct Implementation;
54 
64  template <int structdim, int dim> struct PresentLevelType
65  {
66  struct type
67  {
71  type ()
72  {}
73 
79  type (const int level)
80  {
81  Assert (level == 0, ExcInternalError());
82  (void)level; // removes -Wunused-parameter warning in optimized mode
83  }
84 
90  operator int () const
91  {
92  return 0;
93  }
94 
95  void operator ++ () const
96  {
97  Assert (false, ExcInternalError());
98  }
99 
100  void operator -- () const
101  {
102  Assert (false, ExcInternalError());
103  }
104  };
105  };
106 
107 
117  template <int dim> struct PresentLevelType<dim,dim>
118  {
119  typedef int type;
120  };
121 
122  }
123 }
124 template <int structdim, int dim, int spacedim> class TriaAccessor;
125 template <int dim, int spacedim> class TriaAccessor<0, dim, spacedim>;
126 template <int spacedim> class TriaAccessor<0, 1, spacedim>;
127 
128 // note: the file tria_accessor.templates.h is included at the end of
129 // this file. this includes a lot of templates. originally, this was
130 // only done in debug mode, but led to cyclic reduction problems and
131 // so is now on by default.
132 
133 
138 {
139 //TODO: Write documentation!
143  DeclException0 (ExcCellNotUsed);
155  DeclException0 (ExcCellNotActive);
162  DeclException0 (ExcCellHasNoChildren);
170  DeclException0 (ExcCellHasNoParent);
171 //TODO: Write documentation!
175  DeclException0 (ExcUnusedCellAsChild);
176 //TODO: Write documentation!
180  DeclException1 (ExcCantSetChildren,
181  int,
182  << "You can only set the child index if the cell has no "
183  << "children, or clear it. The given "
184  << "index was " << arg1 << " (-1 means: clear children)");
185 //TODO: Write documentation!
189  DeclException0 (ExcUnusedCellAsNeighbor);
190 //TODO: Write documentation!
194  DeclException0 (ExcUncaughtCase);
195 //TODO: Write documentation!
199  DeclException0 (ExcDereferenceInvalidObject);
200 //TODO: Write documentation!
204  DeclException0 (ExcCantCompareIterators);
205 //TODO: Write documentation!
209  DeclException0 (ExcNeighborIsCoarser);
210 //TODO: Write documentation!
214  DeclException0 (ExcNeighborIsNotCoarser);
226  DeclException0 (ExcFacesHaveNoLevel);
227 //TODO: Write documentation!
231  DeclException1 (ExcSetOnlyEvenChildren,
232  int,
233  << "You can only set the child index of an even numbered child."
234  << "The number of the child given was " << arg1 << ".");
235 }
236 
237 
263 template <int structdim, int dim, int spacedim=dim>
264 class TriaAccessorBase
265 {
266 public:
278  static const unsigned int space_dimension = spacedim;
279 
289  static const unsigned int dimension = dim;
290 
301  static const unsigned int structure_dimension = structdim;
302 
303 protected:
313  typedef void AccessorData;
314 
321  const int level = -1,
322  const int index = -1,
323  const AccessorData * = 0);
324 
330 
342  void copy_from (const TriaAccessorBase &);
343 
349 
356  bool operator < (const TriaAccessorBase &other) const;
357 
358 protected:
376  void operator = (const TriaAccessorBase *);
377 
381  bool operator == (const TriaAccessorBase &) const;
382 
386  bool operator != (const TriaAccessorBase &) const;
387 
406  void operator ++ ();
407 
421  void operator -- ();
432  objects () const;
433 
434 public:
442  typedef void *LocalData;
443 
456  int level () const;
457 
487  int index () const;
488 
498 
506 
510 protected:
516  typename ::internal::TriaAccessor::PresentLevelType<structdim,dim>::type present_level;
517 
525 
531 
532 private:
533 
534  template <typename Accessor> friend class TriaRawIterator;
535  template <typename Accessor> friend class TriaIterator;
536  template <typename Accessor> friend class TriaActiveIterator;
537 };
538 
539 
540 
563 template <int structdim, int dim, int spacedim=dim>
564 class InvalidAccessor : public TriaAccessorBase<structdim,dim,spacedim>
565 {
566 public:
572 
586  InvalidAccessor (const Triangulation<dim,spacedim> *parent = 0,
587  const int level = -1,
588  const int index = -1,
589  const AccessorData *local_data = 0);
590 
605 
613  template <typename OtherAccessor>
614  InvalidAccessor (const OtherAccessor &);
615 
619  void copy_from (const InvalidAccessor &);
620 
624  bool operator == (const InvalidAccessor &) const;
625  bool operator != (const InvalidAccessor &) const;
626 
631  void operator ++ () const;
632  void operator -- () const;
633 
639  bool used () const;
640 
646  bool has_children () const;
647 };
648 
649 
650 
667 template <int structdim, int dim, int spacedim>
668 class TriaAccessor : public TriaAccessorBase<structdim, dim, spacedim>
669 {
670 public:
675  typedef
678 
682  TriaAccessor (const Triangulation<dim,spacedim> *parent = 0,
683  const int level = -1,
684  const int index = -1,
685  const AccessorData *local_data = 0);
686 
708  template <int structdim2, int dim2, int spacedim2>
710 
717  template <int structdim2, int dim2, int spacedim2>
719 
733  bool used () const;
734 
744  int parent_index () const;
745 
746 
777  unsigned int vertex_index (const unsigned int i) const;
778 
783  Point<spacedim> &vertex (const unsigned int i) const;
784 
789  typename ::internal::Triangulation::Iterators<dim,spacedim>::line_iterator
790  line (const unsigned int i) const;
791 
801  unsigned int line_index (const unsigned int i) const;
802 
807  typename ::internal::Triangulation::Iterators<dim,spacedim>::quad_iterator
808  quad (const unsigned int i) const;
809 
819  unsigned int quad_index (const unsigned int i) const;
853  bool face_orientation (const unsigned int face) const;
854 
869  bool face_flip (const unsigned int face) const;
870 
885  bool face_rotation (const unsigned int face) const;
886 
904  bool line_orientation (const unsigned int line) const;
920  bool has_children () const;
921 
928  unsigned int n_children() const;
929 
954  unsigned int number_of_children () const;
955 
977  unsigned int max_refinement_depth () const;
978 
984  child (const unsigned int i) const;
985 
1001  isotropic_child (const unsigned int i) const;
1002 
1008 
1019  int child_index (const unsigned int i) const;
1020 
1028  int isotropic_child_index (const unsigned int i) const;
1052 
1100  void set_boundary_indicator (const types::boundary_id) const;
1101 
1139 
1151  bool at_boundary () const;
1152 
1163  const Boundary<dim,spacedim> &get_boundary () const;
1164 
1180  bool user_flag_set () const;
1181 
1186  void set_user_flag () const;
1187 
1192  void clear_user_flag () const;
1193 
1199  void recursively_set_user_flag () const;
1200 
1206  void recursively_clear_user_flag () const;
1207 
1213  void clear_user_data () const;
1214 
1229  void set_user_pointer (void *p) const;
1230 
1236  void clear_user_pointer () const;
1237 
1260  void *user_pointer () const;
1261 
1299  void recursively_set_user_pointer (void *p) const;
1300 
1310  void recursively_clear_user_pointer () const;
1311 
1325  void set_user_index (const unsigned int p) const;
1326 
1331  void clear_user_index () const;
1332 
1347  unsigned int user_index () const;
1348 
1377  void recursively_set_user_index (const unsigned int p) const;
1378 
1389  void recursively_clear_user_index () const;
1411  double diameter () const;
1412 
1428  double extent_in_direction (const unsigned int axis) const;
1429 
1434  double minimum_vertex_distance () const;
1435 
1448  Point<spacedim> center () const;
1449 
1453  Point<spacedim> barycenter () const;
1454 
1463  double measure () const;
1464 
1491  bool
1493 
1499 private:
1506  void set (const ::internal::Triangulation::TriaObject<structdim> &o) const;
1507 
1518  void set_line_orientation (const unsigned int line,
1519  const bool orientation) const;
1520 
1540  void set_face_orientation (const unsigned int face,
1541  const bool orientation) const;
1542 
1553  void set_face_flip (const unsigned int face,
1554  const bool flip) const;
1555 
1566  void set_face_rotation (const unsigned int face,
1567  const bool rotation) const;
1568 
1574  void set_used_flag () const;
1575 
1581  void clear_used_flag () const;
1582 
1598  void set_refinement_case (const RefinementCase<structdim> &ref_case) const;
1599 
1611  void clear_refinement_case () const;
1612 
1616  void set_parent (const unsigned int parent_index);
1617 
1631  void set_children (const unsigned int i, const int index) const;
1632 
1639  void clear_children () const;
1640 
1641 private:
1659  void operator = (const TriaAccessor &);
1660 
1661  template <int, int> friend class Triangulation;
1662 
1663  friend struct ::internal::Triangulation::Implementation;
1664  friend struct ::internal::TriaAccessor::Implementation;
1665 };
1666 
1667 
1668 
1669 
1670 
1671 
1678 template<int dim, int spacedim>
1679 class TriaAccessor<0, dim, spacedim>
1680 {
1681 private:
1687  TriaAccessor ();
1688 };
1689 
1690 
1691 
1701 template <int spacedim>
1702 class TriaAccessor<0, 1, spacedim>
1703 {
1704 public:
1716  static const unsigned int space_dimension = spacedim;
1717 
1727  static const unsigned int dimension = 1;
1728 
1739  static const unsigned int structure_dimension = 0;
1740 
1744  typedef void AccessorData;
1745 
1753  {
1754  left_vertex,
1755  interior_vertex,
1756  right_vertex
1757  };
1758 
1780  const VertexKind vertex_kind,
1781  const unsigned int vertex_index);
1782 
1793  const int = 0,
1794  const int = 0,
1795  const AccessorData * = 0);
1796 
1802  template <int structdim2, int dim2, int spacedim2>
1804 
1810  template <int structdim2, int dim2, int spacedim2>
1812 
1824  void copy_from (const TriaAccessor &);
1825 
1836 
1842  static int level ();
1843 
1849  int index () const;
1850 
1865  void operator ++ () const;
1866 
1875  void operator -- () const;
1879  bool operator == (const TriaAccessor &) const;
1880 
1884  bool operator != (const TriaAccessor &) const;
1885 
1894  static int parent_index ();
1895 
1927  unsigned int vertex_index (const unsigned int i = 0) const;
1928 
1937  Point<spacedim> &vertex (const unsigned int i = 0) const;
1938 
1945  Point<spacedim> center () const;
1946 
1952  typename ::internal::Triangulation::Iterators<1,spacedim>::line_iterator
1953  static line (const unsigned int);
1954 
1964  static unsigned int line_index (const unsigned int i);
1965 
1970  static
1971  typename ::internal::Triangulation::Iterators<1,spacedim>::quad_iterator
1972  quad (const unsigned int i);
1973 
1983  static unsigned int quad_index (const unsigned int i);
1984 
1996  bool at_boundary () const;
1997 
2015 
2026  static bool face_orientation (const unsigned int face);
2027 
2031  static bool face_flip (const unsigned int face);
2032 
2036  static bool face_rotation (const unsigned int face);
2037 
2041  static bool line_orientation (const unsigned int line);
2042 
2058  static bool has_children ();
2059 
2065  static unsigned int n_children();
2066 
2072  static unsigned int number_of_children ();
2073 
2079  static unsigned int max_refinement_depth ();
2080 
2084  static
2086  child (const unsigned int);
2087 
2091  static
2093  isotropic_child (const unsigned int);
2094 
2098  static
2100 
2104  static
2105  int child_index (const unsigned int i);
2106 
2110  static
2111  int isotropic_child_index (const unsigned int i);
2153  void
2155 
2166  void
2176  bool used () const;
2177 
2178 protected:
2184 
2193 
2199  unsigned int global_vertex_index;
2200 };
2201 
2202 
2203 
2204 
2221 template <int dim, int spacedim=dim>
2222 class CellAccessor : public TriaAccessor<dim,dim,spacedim>
2223 {
2224 public:
2230 
2236 
2248  const int level = -1,
2249  const int index = -1,
2250  const AccessorData *local_data = 0);
2251 
2255  CellAccessor (const TriaAccessor<dim,dim,spacedim> &cell_accessor);
2256 
2278  template <int structdim2, int dim2, int spacedim2>
2280 
2287  template <int structdim2, int dim2, int spacedim2>
2289 
2308  child (const unsigned int i) const;
2309 
2314  TriaIterator<TriaAccessor<dim-1,dim,spacedim> >
2315  face (const unsigned int i) const;
2316 
2328  unsigned int
2329  face_index (const unsigned int i) const;
2330 
2423  neighbor_child_on_subface (const unsigned int face_no,
2424  const unsigned int subface_no) const;
2425 
2438  neighbor (const unsigned int i) const;
2439 
2446  int neighbor_index (const unsigned int i) const;
2447 
2454  int neighbor_level (const unsigned int i) const;
2455 
2477  unsigned int neighbor_of_neighbor (const unsigned int neighbor) const;
2478 
2500  bool neighbor_is_coarser (const unsigned int neighbor) const;
2501 
2520  std::pair<unsigned int, unsigned int>
2521  neighbor_of_coarser_neighbor (const unsigned int neighbor) const;
2522 
2532  unsigned int neighbor_face_no (const unsigned int neighbor) const;
2533 
2553  bool at_boundary (const unsigned int i) const;
2554 
2573  bool at_boundary () const;
2574 
2588  bool has_boundary_lines () const;
2607 
2615 
2619  void clear_refine_flag () const;
2620 
2632  bool flag_for_face_refinement (const unsigned int face_no,
2633  const RefinementCase<dim-1> &face_refinement_case=RefinementCase<dim-1>::isotropic_refinement) const;
2634 
2643  bool flag_for_line_refinement (const unsigned int line_no) const;
2644 
2659  ::internal::SubfaceCase<dim> subface_case(const unsigned int face_no) const;
2660 
2665  bool coarsen_flag_set () const;
2666 
2673  void set_coarsen_flag () const;
2674 
2678  void clear_coarsen_flag () const;
2703 
2716  void set_material_id (const types::material_id new_material_id) const;
2717 
2728  void recursively_set_material_id (const types::material_id new_material_id) const;
2757 
2778  void set_subdomain_id (const types::subdomain_id new_subdomain_id) const;
2779 
2784 
2788  void set_level_subdomain_id (const types::subdomain_id new_level_subdomain_id) const;
2789 
2790 
2809  void recursively_set_subdomain_id (const types::subdomain_id new_subdomain_id) const;
2828  bool direction_flag () const;
2829 
2830 
2831 
2838  parent () const;
2839 
2859  bool active () const;
2860 
2868  bool operator < (const CellAccessor<dim, spacedim> &other) const;
2869 
2870 
2892  bool is_locally_owned () const;
2893 
2898  bool is_locally_owned_on_level () const;
2899 
2933  bool is_ghost () const;
2934 
2977  bool is_artificial () const;
2978 
3002  bool point_inside (const Point<spacedim> &p) const;
3003 
3017  void set_neighbor (const unsigned int i,
3018  const TriaIterator<CellAccessor<dim, spacedim> > &pointer) const;
3019 
3027  CellId id() const;
3028 
3037  DeclException0 (ExcRefineCellNotActive);
3041  DeclException0 (ExcCellFlaggedForRefinement);
3045  DeclException0 (ExcCellFlaggedForCoarsening);
3046 
3047 protected:
3081  unsigned int neighbor_of_neighbor_internal (const unsigned int neighbor) const;
3082 
3089  template<int dim_,int spacedim_ >
3090  bool point_inside_codim(const Point<spacedim_> &p) const;
3091 
3092 
3093 
3094 private:
3102  void set_direction_flag (const bool new_direction_flag) const;
3124 
3125  template <int, int> friend class Triangulation;
3126 
3127  friend struct ::internal::Triangulation::Implementation;
3128 };
3129 
3130 
3131 
3132 /* -------------- declaration of explicit
3133  specializations and general templates ------------- */
3134 
3135 
3136 template <int structdim, int dim, int spacedim>
3137 template <typename OtherAccessor>
3139 InvalidAccessor (const OtherAccessor &)
3140 {
3141  Assert (false,
3142  ExcMessage ("You are attempting an illegal conversion between "
3143  "iterator/accessor types. The constructor you call "
3144  "only exists to make certain template constructs "
3145  "easier to write as dimension independent code but "
3146  "the conversion is not valid in the current context."));
3147 }
3148 
3149 
3150 
3151 template <int structdim, int dim, int spacedim>
3152 template <int structdim2, int dim2, int spacedim2>
3155 {
3156  Assert (false,
3157  ExcMessage ("You are attempting an illegal conversion between "
3158  "iterator/accessor types. The constructor you call "
3159  "only exists to make certain template constructs "
3160  "easier to write as dimension independent code but "
3161  "the conversion is not valid in the current context."));
3162 }
3163 
3164 
3165 
3166 template <int dim, int spacedim>
3167 template <int structdim2, int dim2, int spacedim2>
3170 {
3171  Assert (false,
3172  ExcMessage ("You are attempting an illegal conversion between "
3173  "iterator/accessor types. The constructor you call "
3174  "only exists to make certain template constructs "
3175  "easier to write as dimension independent code but "
3176  "the conversion is not valid in the current context."));
3177 }
3178 
3179 
3180 
3181 template <int structdim, int dim, int spacedim>
3182 template <int structdim2, int dim2, int spacedim2>
3185 {
3186  Assert (false,
3187  ExcMessage ("You are attempting an illegal conversion between "
3188  "iterator/accessor types. The constructor you call "
3189  "only exists to make certain template constructs "
3190  "easier to write as dimension independent code but "
3191  "the conversion is not valid in the current context."));
3192 }
3193 
3194 
3195 
3196 template <int dim, int spacedim>
3197 template <int structdim2, int dim2, int spacedim2>
3200 {
3201  Assert (false,
3202  ExcMessage ("You are attempting an illegal conversion between "
3203  "iterator/accessor types. The constructor you call "
3204  "only exists to make certain template constructs "
3205  "easier to write as dimension independent code but "
3206  "the conversion is not valid in the current context."));
3207 }
3208 
3209 template <int dim, int spacedim>
3210 CellId
3212 {
3213  std::vector<unsigned char> id(this->level(), -1);
3214  unsigned int coarse_index;
3215 
3216  CellAccessor<dim,spacedim> ptr = *this;
3217  while (ptr.level()>0)
3218  {
3219  // find the 'v'st child of our parent we are
3220  unsigned char v=-1;
3221  for (unsigned int c=0; c<ptr.parent()->n_children(); ++c)
3222  {
3223  if (ptr.parent()->child_index(c)==ptr.index())
3224  {
3225  v = c;
3226  break;
3227  }
3228  }
3229 
3230  Assert(v != (unsigned char)-1, ExcInternalError());
3231  id[ptr.level()-1] = v;
3232 
3233  ptr.copy_from( *(ptr.parent()));
3234  }
3235 
3236  Assert(ptr.level()==0, ExcInternalError());
3237  coarse_index = ptr.index();
3238 
3239  return CellId(coarse_index, id);
3240 }
3241 
3242 
3243 #ifndef DOXYGEN
3244 
3245 template <> bool CellAccessor<1,1>::point_inside (const Point<1> &) const;
3246 template <> bool CellAccessor<2,2>::point_inside (const Point<2> &) const;
3247 template <> bool CellAccessor<3,3>::point_inside (const Point<3> &) const;
3248 template <> bool CellAccessor<1,2>::point_inside (const Point<2> &) const;
3249 template <> bool CellAccessor<1,3>::point_inside (const Point<3> &) const;
3250 template <> bool CellAccessor<2,3>::point_inside (const Point<3> &) const;
3251 // -------------------------------------------------------------------
3252 
3253 #endif // DOXYGEN
3254 
3255 DEAL_II_NAMESPACE_CLOSE
3256 
3257 // include more templates in debug and optimized mode
3258 # include "tria_accessor.templates.h"
3259 
3260 #endif
void set_subdomain_id(const types::subdomain_id new_subdomain_id) const
bool flag_for_line_refinement(const unsigned int line_no) const
void set_neighbor(const unsigned int i, const TriaIterator< CellAccessor< dim, spacedim > > &pointer) const
const Triangulation< dim, spacedim > & get_triangulation() const
TriaAccessor< dim, dim, spacedim >::AccessorData AccessorData
TriaIterator< CellAccessor< dim, spacedim > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
unsigned int user_index() const
void set_user_index(const unsigned int p) const
void set_user_pointer(void *p) const
void set(const ::internal::Triangulation::TriaObject< structdim > &o) const
void clear_children() const
void clear_user_flag() const
bool operator<(const TriaAccessorBase &other) const
TriaIterator< CellAccessor< dim, spacedim > > child(const unsigned int i) const
void clear_user_pointer() const
std::pair< unsigned int, unsigned int > neighbor_of_coarser_neighbor(const unsigned int neighbor) const
unsigned char material_id
Definition: types.h:142
void recursively_set_user_pointer(void *p) const
double measure() const
TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
TriaAccessorBase & operator=(const TriaAccessorBase &)
::ExceptionBase & ExcMessage(std::string arg1)
unsigned int neighbor_of_neighbor_internal(const unsigned int neighbor) const
CellId id() const
bool face_orientation(const unsigned int face) const
unsigned int neighbor_face_no(const unsigned int neighbor) const
TriaIterator< TriaAccessor< dim-1, dim, spacedim > > face(const unsigned int i) const
unsigned int neighbor_of_neighbor(const unsigned int neighbor) const
void copy_from(const TriaAccessorBase &)
::internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
TriaIterator< CellAccessor< dim, spacedim > > neighbor(const unsigned int i) const
void set_face_flip(const unsigned int face, const bool flip) const
bool line_orientation(const unsigned int line) const
int neighbor_level(const unsigned int i) const
double minimum_vertex_distance() const
const Triangulation< 1, spacedim > * tria
const Boundary< dim, spacedim > & get_boundary() const
typename::internal::Triangulation::Iterators< dim, spacedim >::line_iterator line(const unsigned int i) const
types::subdomain_id level_subdomain_id() const
bool face_flip(const unsigned int face) const
void set_user_flag() const
bool has_children() const
bool is_locally_owned() const
void copy_from(const InvalidAccessor &)
bool is_translation_of(const TriaIterator< TriaAccessor< structdim, dim, spacedim > > &o) const
void recursively_clear_user_index() const
static const unsigned int space_dimension
unsigned int max_refinement_depth() const
void clear_user_data() const
void clear_used_flag() const
bool at_boundary() const
TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
void set_refinement_case(const RefinementCase< structdim > &ref_case) const
double diameter() const
bool point_inside_codim(const Point< spacedim_ > &p) const
void set_used_flag() const
void operator=(const CellAccessor< dim, spacedim > &)
void set_parent(const unsigned int parent_index)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
void set_direction_flag(const bool new_direction_flag) const
#define Assert(cond, exc)
Definition: exceptions.h:299
bool at_boundary() const
types::boundary_id boundary_indicator() const
void set_boundary_indicator(const types::boundary_id) const
TriaAccessor(const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *local_data=0)
unsigned int vertex_index(const unsigned int i) const
bool user_flag_set() const
const Triangulation< dim, spacedim > * tria
#define DeclException0(Exception0)
Definition: exceptions.h:505
types::material_id material_id() const
void recursively_clear_user_flag() const
bool is_locally_owned_on_level() const
::internal::Triangulation::TriaObjects<::internal::Triangulation::TriaObject< structdim > > & objects() const
bool is_artificial() const
bool coarsen_flag_set() const
TriaIterator< TriaAccessor< structdim, dim, spacedim > > isotropic_child(const unsigned int i) const
void set_all_boundary_indicators(const types::boundary_id) const
typename::internal::TriaAccessor::PresentLevelType< structdim, dim >::type present_level
TriaIterator< CellAccessor< dim, spacedim > > parent() const
void clear_refine_flag() const
void set_refine_flag(const RefinementCase< dim > ref_case=RefinementCase< dim >::isotropic_refinement) const
bool operator==(const TriaAccessorBase &) const
unsigned int face_index(const unsigned int i) const
unsigned int subdomain_id
Definition: types.h:43
void set_level_subdomain_id(const types::subdomain_id new_level_subdomain_id) const
void clear_coarsen_flag() const
typename::internal::Triangulation::Iterators< dim, spacedim >::quad_iterator quad(const unsigned int i) const
bool point_inside(const Point< spacedim > &p) const
Point< spacedim > center() const
int neighbor_index(const unsigned int i) const
bool has_boundary_lines() const
void recursively_set_material_id(const types::material_id new_material_id) const
Definition: cell_id.h:40
RefinementCase< dim > refine_flag_set() const
bool neighbor_is_coarser(const unsigned int neighbor) const
void set_face_orientation(const unsigned int face, const bool orientation) const
unsigned int number_of_children() const
TriaIterator< TriaAccessor< structdim, dim, spacedim > > child(const unsigned int i) const
unsigned int quad_index(const unsigned int i) const
void clear_refinement_case() 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 n_children() const
IteratorState::IteratorStates state() const
void recursively_set_user_index(const unsigned int p) const
void set_coarsen_flag() const
static const unsigned int structure_dimension
CellAccessor(const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *local_data=0)
void set_face_rotation(const unsigned int face, const bool rotation) const
int isotropic_child_index(const unsigned int i) const
void set_children(const unsigned int i, const int index) const
Triangulation< dim, spacedim > Container
void * user_pointer() const
void recursively_clear_user_pointer() const
void recursively_set_subdomain_id(const types::subdomain_id new_subdomain_id) const
int child_index(const unsigned int i) const
InvalidAccessor(const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *local_data=0)
double extent_in_direction(const unsigned int axis) const
unsigned char boundary_id
Definition: types.h:128
types::subdomain_id subdomain_id() const
void recursively_set_user_flag() const
static const unsigned int dimension
void clear_user_index() const
TriaAccessorBase(const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *=0)
::ExceptionBase & ExcInternalError()
void set_material_id(const types::material_id new_material_id) const
void operator=(const TriaAccessor &)
bool operator==(const InvalidAccessor &) const
Point< spacedim > barycenter() const
bool face_rotation(const unsigned int face) const
RefinementCase< structdim > refinement_case() const
void set_line_orientation(const unsigned int line, const bool orientation) const
bool direction_flag() const
bool operator!=(const TriaAccessorBase &) const
Point< spacedim > & vertex(const unsigned int i) const
unsigned int line_index(const unsigned int i) const