Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
tria_accessor.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2022 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
31
33#include <boost/container/small_vector.hpp>
35
36#include <utility>
37
38
40
41// Forward declarations
42#ifndef DOXYGEN
43template <int dim, int spacedim>
44class Triangulation;
45template <typename Accessor>
46class TriaRawIterator;
47template <typename Accessor>
48class TriaIterator;
49template <typename Accessor>
51
52namespace parallel
53{
54 template <int dim, int spacedim>
55 class TriangulationBase;
56}
57
58template <int dim, int spacedim>
59class Manifold;
60
61template <int dim, int spacedim>
62class Mapping;
63#endif
64
65namespace 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 {
99 (void)level; // removes -Wunused-parameter warning in optimized mode
100 }
101
105 operator int() const
106 {
107 return 0;
108 }
109
110 void
112 {
113 Assert(false, ExcInternalError());
114 }
115
116 void
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
137template <int structdim, int dim, int spacedim>
138class TriaAccessor;
139template <int dim, int spacedim>
140class TriaAccessor<0, dim, spacedim>;
141template <int spacedim>
142class 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
300template <int structdim, int dim, int spacedim = dim>
302{
303public:
309 static constexpr unsigned int space_dimension = spacedim;
310
316 static constexpr unsigned int dimension = dim;
317
323 static const unsigned int structure_dimension = structdim;
324
334 void
335 operator=(const TriaAccessorBase *) = delete;
336
337protected:
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
387protected:
391 bool
393
397 bool
399
413 void
415
423 void
433 objects() const;
434
435public:
441 using LocalData = void *;
442
466 int
467 level() const;
468
495 int
496 index() const;
497
503 state() const;
504
511
515protected:
520 typename ::internal::TriaAccessorImplementation::
521 PresentLevelType<structdim, dim>::type present_level;
522
528
533
534private:
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
565template <int structdim, int dim, int spacedim = dim>
566class InvalidAccessor : public TriaAccessorBase<structdim, dim, spacedim>
567{
568public:
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
614 bool
616
620 void
621 operator++() const;
622 void
623 operator--() const;
624
629 bool
630 used() const;
631
636 bool
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
705template <int structdim, int dim, int spacedim>
706class TriaAccessor : public TriaAccessorBase<structdim, dim, spacedim>
707{
708public:
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
766 operator=(const TriaAccessor &) = delete;
767
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;
906 unsigned char
907 combined_face_orientation(const unsigned int face) const;
908
920 bool
921 face_orientation(const unsigned int face) const;
922
932 bool
933 face_flip(const unsigned int face) const;
934
944 bool
945 face_rotation(const unsigned int face) const;
946
959 bool
960 line_orientation(const unsigned int line) const;
975 bool
977
982 unsigned int
983 n_children() const;
984
989 unsigned int
991
1005 unsigned int
1007
1021 unsigned int
1023
1028 child(const unsigned int i) const;
1029
1034 unsigned int
1037
1047 isotropic_child(const unsigned int i) const;
1048
1054
1060 int
1061 child_index(const unsigned int i) const;
1062
1068 int
1069 isotropic_child_index(const unsigned int i) const;
1093
1123 void
1125
1156 void
1158
1166 bool
1168
1180
1203
1221 void
1223
1237 void
1239
1256 bool
1258
1264 void
1266
1272 void
1274
1280 void
1282
1288 void
1290
1296 void
1298
1310 void
1311 set_user_pointer(void *p) const;
1312
1318 void
1320
1336 void *
1338
1360 void
1362
1369 void
1371
1381 void
1382 set_user_index(const unsigned int p) const;
1383
1389 void
1391
1403 unsigned int
1404 user_index() const;
1405
1423 void
1424 recursively_set_user_index(const unsigned int p) const;
1425
1434 void
1470 double
1471 diameter() const;
1472
1499 std::pair<Point<spacedim>, double>
1501
1511 bounding_box() const;
1512
1522 double
1523 extent_in_direction(const unsigned int axis) const;
1524
1528 double
1530
1545 intermediate_point(const Point<structdim> &coordinates) const;
1546
1571
1607 center(const bool respect_manifold = false,
1608 const bool interpolate_from_surrounding = false) const;
1609
1628 barycenter() const;
1629
1655 double
1656 measure() const;
1657
1672 bool
1675
1681
1685 unsigned int
1686 n_vertices() const;
1687
1691 unsigned int
1692 n_lines() const;
1693
1699 unsigned int
1700 n_faces() const;
1701
1708
1715
1724
1729private:
1734 void
1736
1744 void
1746 const std::initializer_list<int> &new_indices) const;
1747
1751 void
1753 const std::initializer_list<unsigned int> &new_indices) const;
1754
1762 void
1763 set_line_orientation(const unsigned int line, const bool orientation) const;
1764
1775 void
1776 set_face_orientation(const unsigned int face, const bool orientation) const;
1777
1784 void
1785 set_face_flip(const unsigned int face, const bool flip) const;
1786
1793 void
1794 set_face_rotation(const unsigned int face, const bool rotation) const;
1795
1799 void
1801
1805 void
1807
1816 void
1818
1826 void
1828
1835 void
1836 set_children(const unsigned int i, const int index) const;
1837
1842 void
1844
1845private:
1846 template <int, int>
1847 friend class Triangulation;
1848
1849 friend struct ::internal::TriangulationImplementation::Implementation;
1850 friend struct ::internal::TriangulationImplementation::
1851 ImplementationMixedMesh;
1852 friend struct ::internal::TriaAccessorImplementation::Implementation;
1853};
1854
1855
1856
1875template <int dim, int spacedim>
1876class TriaAccessor<0, dim, spacedim>
1877{
1878public:
1884 static constexpr unsigned int space_dimension = spacedim;
1885
1891 static constexpr unsigned int dimension = dim;
1892
1898 static const unsigned int structure_dimension = 0;
1899
1903 using AccessorData = void;
1904
1910 const unsigned int vertex_index);
1911
1918 const int level = 0,
1919 const int index = 0,
1920 const AccessorData * = nullptr);
1921
1925 template <int structdim2, int dim2, int spacedim2>
1927
1931 template <int structdim2, int dim2, int spacedim2>
1933
1938 state() const;
1939
1944 static int
1946
1951 int
1952 index() const;
1953
1960
1970 void
1972
1976 void
1981 bool
1982 operator==(const TriaAccessor &) const;
1983
1987 bool
1988 operator!=(const TriaAccessor &) const;
1989
2017 unsigned int
2018 vertex_index(const unsigned int i = 0) const;
2019
2026 vertex(const unsigned int i = 0) const;
2027
2032 typename ::internal::TriangulationImplementation::
2033 Iterators<dim, spacedim>::line_iterator static line(const unsigned int);
2034
2038 static unsigned int
2039 line_index(const unsigned int i);
2040
2044 static typename ::internal::TriangulationImplementation::
2045 Iterators<dim, spacedim>::quad_iterator
2046 quad(const unsigned int i);
2047
2051 static unsigned int
2052 quad_index(const unsigned int i);
2053
2069 double
2070 diameter() const;
2071
2079 double
2080 extent_in_direction(const unsigned int axis) const;
2081
2090 center(const bool respect_manifold = false,
2091 const bool interpolate_from_surrounding = false) const;
2092
2101 double
2102 measure() const;
2117 static unsigned char
2118 combined_face_orientation(const unsigned int face);
2119
2123 static bool
2124 face_orientation(const unsigned int face);
2125
2129 static bool
2130 face_flip(const unsigned int face);
2131
2135 static bool
2136 face_rotation(const unsigned int face);
2137
2141 static bool
2142 line_orientation(const unsigned int line);
2143
2158 static bool
2160
2165 static unsigned int
2167
2172 static unsigned int
2174
2179 static unsigned int
2181
2182
2186 static unsigned int
2188
2192 static unsigned int
2194
2199 child(const unsigned int);
2200
2205 isotropic_child(const unsigned int);
2206
2210 static RefinementCase<0>
2212
2216 static int
2217 child_index(const unsigned int i);
2218
2222 static int
2223 isotropic_child_index(const unsigned int i);
2231 bool
2232 used() const;
2233
2234protected:
2242 void
2244
2253 bool
2254 operator<(const TriaAccessor &other) const;
2255
2260
2265
2266private:
2267 template <typename Accessor>
2268 friend class TriaRawIterator;
2269 template <typename Accessor>
2270 friend class TriaIterator;
2271 template <typename Accessor>
2273};
2274
2275
2276
2293template <int spacedim>
2294class TriaAccessor<0, 1, spacedim>
2295{
2296public:
2302 static constexpr unsigned int space_dimension = spacedim;
2303
2309 static constexpr unsigned int dimension = 1;
2310
2316 static const unsigned int structure_dimension = 0;
2317
2321 using AccessorData = void;
2322
2328 {
2340 right_vertex
2342
2355 const VertexKind vertex_kind,
2356 const unsigned int vertex_index);
2357
2364 const int = 0,
2365 const int = 0,
2366 const AccessorData * = nullptr);
2367
2371 template <int structdim2, int dim2, int spacedim2>
2373
2377 template <int structdim2, int dim2, int spacedim2>
2379
2384 void
2386
2392 void
2394
2402
2407 static int
2409
2414 int
2415 index() const;
2416
2423
2434 void
2435 operator++() const;
2436
2441 void
2442 operator--() const;
2446 bool
2447 operator==(const TriaAccessor &) const;
2448
2452 bool
2453 operator!=(const TriaAccessor &) const;
2454
2463 bool
2464 operator<(const TriaAccessor &other) const;
2465
2492 unsigned int
2493 vertex_index(const unsigned int i = 0) const;
2494
2501 vertex(const unsigned int i = 0) const;
2502
2508 center() const;
2509
2514 typename ::internal::TriangulationImplementation::
2515 Iterators<1, spacedim>::line_iterator static line(const unsigned int);
2516
2523 static unsigned int
2524 line_index(const unsigned int i);
2525
2529 static typename ::internal::TriangulationImplementation::
2530 Iterators<1, spacedim>::quad_iterator
2531 quad(const unsigned int i);
2532
2539 static unsigned int
2540 quad_index(const unsigned int i);
2541
2551 bool
2553
2570
2574 const Manifold<1, spacedim> &
2576
2585
2586
2598 bool
2599 user_flag_set() const;
2600
2606 void
2607 set_user_flag() const;
2608
2614 void
2615 clear_user_flag() const;
2616
2622 void
2624
2630 void
2632
2638 void
2639 clear_user_data() const;
2640
2652 void
2653 set_user_pointer(void *p) const;
2654
2660 void
2661 clear_user_pointer() const;
2662
2678 void *
2679 user_pointer() const;
2680
2702 void
2703 recursively_set_user_pointer(void *p) const;
2704
2711 void
2713
2723 void
2724 set_user_index(const unsigned int p) const;
2725
2731 void
2732 clear_user_index() const;
2733
2745 unsigned int
2746 user_index() const;
2747
2765 void
2766 recursively_set_user_index(const unsigned int p) const;
2767
2776 void
2792 static unsigned char
2793 combined_face_orientation(const unsigned int face);
2794
2798 static bool
2799 face_orientation(const unsigned int face);
2800
2804 static bool
2805 face_flip(const unsigned int face);
2806
2810 static bool
2811 face_rotation(const unsigned int face);
2812
2816 static bool
2817 line_orientation(const unsigned int line);
2818
2833 static bool
2835
2840 static unsigned int
2842
2847 static unsigned int
2849
2854 static unsigned int
2856
2857
2861 static unsigned int
2863
2867 static unsigned int
2869
2874 child(const unsigned int);
2875
2880 isotropic_child(const unsigned int);
2881
2885 static RefinementCase<0>
2887
2891 static int
2892 child_index(const unsigned int i);
2893
2897 static int
2898 isotropic_child_index(const unsigned int i);
2931 void
2933
2940 void
2942
2954 void
2956
2968 void
2977 bool
2978 used() const;
2979
2985
2989 unsigned int
2990 n_vertices() const;
2991
2995 unsigned int
2996 n_lines() const;
2997
3004
3011
3012protected:
3017
3023
3028};
3029
3030
3031
3047template <int dim, int spacedim = dim>
3048class CellAccessor : public TriaAccessor<dim, dim, spacedim>
3049{
3050public:
3055
3060
3072 const int level = -1,
3073 const int index = -1,
3074 const AccessorData * local_data = nullptr);
3075
3080
3093 template <int structdim2, int dim2, int spacedim2>
3095
3100 template <int structdim2, int dim2, int spacedim2>
3102
3107
3111 // NOLINTNEXTLINE OSX does not compile with noexcept
3113
3117 ~CellAccessor() = default;
3118
3130
3134 // NOLINTNEXTLINE OSX does not compile with noexcept
3137
3154 child(const unsigned int i) const;
3155
3159 boost::container::small_vector<TriaIterator<CellAccessor<dim, spacedim>>,
3162
3166 TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>
3167 face(const unsigned int i) const;
3168
3173 unsigned int
3176
3180 boost::container::small_vector<
3181 TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>,
3184
3194 unsigned int
3195 face_index(const unsigned int i) const;
3196
3245 neighbor_child_on_subface(const unsigned int face_no,
3246 const unsigned int subface_no) const;
3247
3273 neighbor(const unsigned int face_no) const;
3274
3282 int
3283 neighbor_index(const unsigned int face_no) const;
3284
3292 int
3293 neighbor_level(const unsigned int face_no) const;
3294
3306 unsigned int
3307 neighbor_of_neighbor(const unsigned int face_no) const;
3308
3319 bool
3320 neighbor_is_coarser(const unsigned int face_no) const;
3321
3336 std::pair<unsigned int, unsigned int>
3337 neighbor_of_coarser_neighbor(const unsigned int neighbor) const;
3338
3345 unsigned int
3346 neighbor_face_no(const unsigned int neighbor) const;
3347
3351 static bool
3353
3367 bool
3368 has_periodic_neighbor(const unsigned int i) const;
3369
3387 periodic_neighbor(const unsigned int i) const;
3388
3397 neighbor_or_periodic_neighbor(const unsigned int i) const;
3398
3414 periodic_neighbor_child_on_subface(const unsigned int face_no,
3415 const unsigned int subface_no) const;
3416
3427 std::pair<unsigned int, unsigned int>
3428 periodic_neighbor_of_coarser_periodic_neighbor(const unsigned face_no) const;
3429
3435 int
3436 periodic_neighbor_index(const unsigned int i) const;
3437
3443 int
3444 periodic_neighbor_level(const unsigned int i) const;
3445
3460 unsigned int
3461 periodic_neighbor_of_periodic_neighbor(const unsigned int i) const;
3462
3468 unsigned int
3469 periodic_neighbor_face_no(const unsigned int i) const;
3470
3477 bool
3478 periodic_neighbor_is_coarser(const unsigned int i) const;
3479
3496 bool
3497 at_boundary(const unsigned int i) const;
3498
3507 bool
3508 at_boundary() const;
3509
3517 bool
3518 has_boundary_lines() const;
3546
3564 void
3567
3571 void
3573
3581 bool
3583 const unsigned int face_no,
3584 const RefinementCase<dim - 1> &face_refinement_case =
3586
3592 bool
3593 flag_for_line_refinement(const unsigned int line_no) const;
3594
3604 subface_case(const unsigned int face_no) const;
3605
3609 bool
3611
3616 void
3618
3622 void
3647 material_id() const;
3648
3660 void
3661 set_material_id(const types::material_id new_material_id) const;
3662
3671 void
3672 recursively_set_material_id(const types::material_id new_material_id) const;
3700
3716 void
3717 set_subdomain_id(const types::subdomain_id new_subdomain_id) const;
3718
3729
3734 void
3736 const types::subdomain_id new_level_subdomain_id) const;
3737
3738
3754 void
3756 const types::subdomain_id new_subdomain_id) const;
3780
3788
3802 bool
3803 direction_flag() const;
3804
3830 unsigned int
3832
3840 int
3841 parent_index() const;
3842
3849 parent() const;
3850
3870 bool
3871 is_active() const;
3872
3892 bool
3894
3899 bool
3901
3935 bool
3936 is_ghost() const;
3937
3943 bool
3945
3972 bool
3974
3981 bool
3983
3997 bool
3999
4008 void
4009 set_neighbor(const unsigned int i,
4010 const TriaIterator<CellAccessor<dim, spacedim>> &pointer) const;
4011
4025 CellId
4026 id() const;
4027
4028 using TriaAccessor<dim, dim, spacedim>::diameter;
4029
4033 double
4034 diameter(const Mapping<dim, spacedim> &mapping) const;
4035
4053
4054protected:
4070 unsigned int
4071 neighbor_of_neighbor_internal(const unsigned int neighbor) const;
4072
4078 template <int dim_, int spacedim_>
4079 bool
4080 point_inside_codim(const Point<spacedim_> &p) const;
4081
4082
4083
4084private:
4089 void
4090 set_active_cell_index(const unsigned int active_cell_index) const;
4091
4095 void
4097
4101 void
4103
4107 void
4108 set_parent(const unsigned int parent_index);
4109
4116 void
4117 set_direction_flag(const bool new_direction_flag) const;
4118
4119 template <int, int>
4120 friend class Triangulation;
4121
4122 template <int, int>
4124
4125 friend struct ::internal::TriangulationImplementation::Implementation;
4126 friend struct ::internal::TriangulationImplementation::
4127 ImplementationMixedMesh;
4128};
4129
4130
4131
4132/* ----- declaration of explicit specializations and general templates ----- */
4133
4134
4135template <int structdim, int dim, int spacedim>
4136template <typename OtherAccessor>
4138 const OtherAccessor &)
4139{
4140 Assert(false,
4141 ExcMessage("You are attempting an illegal conversion between "
4142 "iterator/accessor types. The constructor you call "
4143 "only exists to make certain template constructs "
4144 "easier to write as dimension independent code but "
4145 "the conversion is not valid in the current context."));
4146}
4147
4148
4149
4150template <int structdim, int dim, int spacedim>
4151template <int structdim2, int dim2, int spacedim2>
4154{
4155 Assert(false,
4156 ExcMessage("You are attempting an illegal conversion between "
4157 "iterator/accessor types. The constructor you call "
4158 "only exists to make certain template constructs "
4159 "easier to write as dimension independent code but "
4160 "the conversion is not valid in the current context."));
4161}
4162
4163
4164
4165template <int dim, int spacedim>
4166template <int structdim2, int dim2, int spacedim2>
4169{
4170 Assert(false,
4171 ExcMessage("You are attempting an illegal conversion between "
4172 "iterator/accessor types. The constructor you call "
4173 "only exists to make certain template constructs "
4174 "easier to write as dimension independent code but "
4175 "the conversion is not valid in the current context."));
4176}
4177
4178
4179
4180template <int structdim, int dim, int spacedim>
4181template <int structdim2, int dim2, int spacedim2>
4184{
4185 Assert(false,
4186 ExcMessage("You are attempting an illegal conversion between "
4187 "iterator/accessor types. The constructor you call "
4188 "only exists to make certain template constructs "
4189 "easier to write as dimension independent code but "
4190 "the conversion is not valid in the current context."));
4191}
4192
4193
4194
4195template <int dim, int spacedim>
4196template <int structdim2, int dim2, int spacedim2>
4199{
4200 Assert(false,
4201 ExcMessage("You are attempting an illegal conversion between "
4202 "iterator/accessor types. The constructor you call "
4203 "only exists to make certain template constructs "
4204 "easier to write as dimension independent code but "
4205 "the conversion is not valid in the current context."));
4206}
4207
4208
4209#ifndef DOXYGEN
4210
4211template <>
4212bool
4214template <>
4215bool
4217template <>
4218bool
4220template <>
4221bool
4223template <>
4224bool
4226template <>
4227bool
4229// -------------------------------------------------------------------
4230
4231template <>
4232void
4234
4235#endif // DOXYGEN
4236
4238
4239// include more templates in debug and optimized mode
4240#include "tria_accessor.templates.h"
4241
4242#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
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_neighbor(const unsigned int i, const TriaIterator< CellAccessor< dim, spacedim > > &pointer) const
TriaIterator< CellAccessor< dim, spacedim > > neighbor(const unsigned int face_no) const
bool is_artificial_on_level() 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
RefinementCase< dim > refine_flag_set() 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
CellAccessor< dim, spacedim > & operator=(CellAccessor< dim, spacedim > &&)=default
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
bool neighbor_is_coarser(const unsigned int face_no) const
TriaIterator< CellAccessor< dim, spacedim > > child(const unsigned int i) const
void set_global_level_cell_index(const types::global_cell_index index) const
bool has_periodic_neighbor(const unsigned int i) 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
CellAccessor< dim, spacedim > & operator=(const CellAccessor< dim, spacedim > &)=delete
boost::container::small_vector< TriaIterator< TriaAccessor< dim - 1, dim, spacedim > >, GeometryInfo< dim >::faces_per_cell > face_iterators() const
unsigned int neighbor_of_neighbor(const unsigned int face_no) const
unsigned int face_iterator_to_index(const TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > &face) const
void set_material_id(const types::material_id new_material_id) const
bool is_locally_owned() const
void set_refine_flag(const RefinementCase< dim > ref_case=RefinementCase< dim >::isotropic_refinement) 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
boost::container::small_vector< TriaIterator< CellAccessor< dim, spacedim > >, GeometryInfo< dim >::max_children_per_cell > child_iterators() 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
CellId id() const
::internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
bool is_artificial() const
bool is_ghost_on_level() 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
Definition: cell_id.h:71
void set_manifold_id(const types::manifold_id) const
InvalidAccessor(const InvalidAccessor &)
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 operator==(const InvalidAccessor &) const
Point< spacedim > & vertex(const unsigned int i) const
void operator++() const
bool operator!=(const InvalidAccessor &) const
bool has_children() const
unsigned int user_index() const
typename::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(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:311
Definition: point.h:111
bool operator!=(const TriaAccessorBase &) const
static constexpr unsigned int space_dimension
static constexpr unsigned int dimension
TriaAccessorBase(const TriaAccessorBase &)
void operator=(const TriaAccessorBase *)=delete
static const unsigned int structure_dimension
TriaAccessorBase(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *=nullptr)
void copy_from(const TriaAccessorBase &)
const Triangulation< dim, spacedim > & get_triangulation() const
IteratorState::IteratorStates state() const
int index() const
bool operator<(const TriaAccessorBase &other) const
typename::internal::TriaAccessorImplementation::PresentLevelType< structdim, dim >::type present_level
::internal::TriangulationImplementation::TriaObjects & objects() const
TriaAccessorBase & operator=(const TriaAccessorBase &)
int level() const
const Triangulation< dim, spacedim > * tria
bool operator==(const TriaAccessorBase &) const
static TriaIterator< TriaAccessor< 0, 1, spacedim > > isotropic_child(const unsigned int)
Return an invalid object.
static unsigned int line_index(const unsigned int i)
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
static bool face_orientation(const unsigned int face)
Always return false.
static unsigned int number_of_children()
static unsigned int child_iterator_to_index(const TriaIterator< TriaAccessor< 0, 1, spacedim > > &)
Return an invalid unsigned integer.
static unsigned int max_refinement_depth()
static unsigned int quad_index(const unsigned int i)
static typename::internal::TriangulationImplementation::Iterators< 1, spacedim >::quad_iterator quad(const unsigned int i)
const Manifold< 1, spacedim > & get_manifold() const
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 > > child(const unsigned int)
Return an invalid object.
unsigned int n_lines() const
const Triangulation< 1, spacedim > & get_triangulation() 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 bool line_orientation(const unsigned int line)
Always return false.
static typename::internal::TriangulationImplementation::Iterators< 1, spacedim >::line_iterator line(const unsigned int)
static RefinementCase< 0 > refinement_case()
ReferenceCell reference_cell() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > line_indices() const
unsigned int vertex_index(const unsigned int i=0) const
void copy_from(const TriaAccessor &)
void set_all_boundary_ids(const types::boundary_id) const
static IteratorState::IteratorStates state()
static bool face_rotation(const unsigned int face)
Always return false.
Point< spacedim > center() const
static int child_index(const unsigned int i)
Returns -1.
types::manifold_id manifold_id() const
static unsigned int n_children()
TriaAccessor(const TriaAccessor< structdim2, dim2, spacedim2 > &)
static unsigned char combined_face_orientation(const unsigned int face)
Always return 0.
void set_boundary_id(const types::boundary_id) const
Point< spacedim > & vertex(const unsigned int i=0) const
unsigned int n_vertices() const
void set_manifold_id(const types::manifold_id)
void copy_from(const TriaAccessorBase< 0, 1, spacedim > &)
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)
Point< spacedim > & vertex(const unsigned int i=0) const
static typename::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(const unsigned int)
bool operator!=(const TriaAccessor &) const
static RefinementCase< 0 > refinement_case()
void copy_from(const TriaAccessor &)
static TriaIterator< TriaAccessor< 0, dim, spacedim > > child(const unsigned int)
Return an invalid object.
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)
static TriaIterator< TriaAccessor< 0, dim, spacedim > > isotropic_child(const unsigned int)
Return an invalid object.
static typename::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator quad(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.
TriaAccessor(const Triangulation< dim, spacedim > *tria=nullptr, const int level=0, const int index=0, const AccessorData *=nullptr)
static unsigned char combined_face_orientation(const unsigned int face)
Always return 0.
static unsigned int n_children()
const Triangulation< dim, spacedim > & get_triangulation() const
static unsigned int child_iterator_to_index(const TriaIterator< TriaAccessor< 0, dim, spacedim > > &)
Return an invalid unsigned integer.
static bool face_flip(const unsigned int face)
Always return false.
static bool line_orientation(const unsigned int line)
Always return false.
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.
static bool face_orientation(const unsigned int face)
Always return false.
static unsigned int n_active_descendants()
TriaAccessor(const TriaAccessor< structdim2, dim2, spacedim2 > &)
Point< spacedim > center(const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
IteratorState::IteratorStates state() const
static unsigned int number_of_children()
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
void clear_children() 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 n_active_descendants() const
unsigned char combined_face_orientation(const unsigned int face) 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
TriaIterator< TriaAccessor< 0, dim, spacedim > > vertex_iterator(const unsigned int i) const
void clear_user_data() const
TriaAccessor(const TriaAccessor< structdim2, dim2, spacedim2 > &)
Point< structdim > real_to_unit_cell_affine_approximation(const Point< spacedim > &point) const
void recursively_clear_user_index() const
unsigned int line_index(const unsigned int i) const
bool face_orientation(const unsigned int face) const
const Manifold< dim, spacedim > & get_manifold() 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
std_cxx20::ranges::iota_view< unsigned int, unsigned int > line_indices() const
typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData AccessorData
Point< spacedim > barycenter() const
BoundingBox< spacedim > bounding_box() const
void clear_user_flag() const
TriaIterator< TriaAccessor< structdim, dim, spacedim > > child(const unsigned int i) const
unsigned int n_children() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices() const
void recursively_set_user_flag() const
void set_boundary_id(const types::boundary_id) const
bool is_translation_of(const TriaIterator< TriaAccessor< structdim, dim, spacedim > > &o) 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
void clear_used_flag() const
std::pair< Point< spacedim >, double > enclosing_ball() const
unsigned int vertex_index(const unsigned int i) const
void * user_pointer() const
int isotropic_child_index(const unsigned int i) 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
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
Point< spacedim > & vertex(const unsigned int i) const
void set_line_orientation(const unsigned int line, const bool orientation) const
unsigned int quad_index(const unsigned int i) const
TriaAccessor & operator=(const TriaAccessor &)=delete
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
RefinementCase< structdim > refinement_case() const
typename::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(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)
unsigned int child_iterator_to_index(const TriaIterator< TriaAccessor< structdim, dim, spacedim > > &child) 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
TriaAccessor & operator=(TriaAccessor &&)=default
TriaIterator< TriaAccessor< structdim, dim, spacedim > > isotropic_child(const unsigned int i) const
typename::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator quad(const unsigned int i) const
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:456
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:495
Point< 3 > center
unsigned int level
Definition: grid_out.cc:4606
static ::ExceptionBase & ExcCantCompareIterators()
#define DeclException0(Exception0)
Definition: exceptions.h:464
static ::ExceptionBase & ExcCellHasNoParent()
static ::ExceptionBase & ExcCellFlaggedForRefinement()
static ::ExceptionBase & ExcCellHasNoChildren()
static ::ExceptionBase & ExcCellNotUsed()
static ::ExceptionBase & ExcNoPeriodicNeighbor()
static ::ExceptionBase & ExcCantSetChildren(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & ExcNeighborIsNotCoarser()
static ::ExceptionBase & ExcSetOnlyEvenChildren(int arg1)
static ::ExceptionBase & ExcNeighborIsCoarser()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcFacesHaveNoLevel()
static ::ExceptionBase & ExcDereferenceInvalidObject(AccessorType arg1)
static ::ExceptionBase & ExcCellFlaggedForCoarsening()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
static ::ExceptionBase & ExcRefineCellNotActive()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcCellNotActive()
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.
boost::integer_range< IncrementableType > iota_view
Definition: iota_view.h:46
unsigned int manifold_id
Definition: types.h:141