deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_tria_accessor_h
16#define dealii_tria_accessor_h
17
18
19#include <deal.II/base/config.h>
20
24#include <deal.II/base/point.h>
25
30
31#include <boost/container/small_vector.hpp>
32
33#include <utility>
34
35
37
38// Forward declarations
39#ifndef DOXYGEN
40template <int dim, int spacedim>
42class Triangulation;
43template <typename Accessor>
44class TriaRawIterator;
45template <typename Accessor>
46class TriaIterator;
47template <typename Accessor>
49
50namespace parallel
51{
52 template <int dim, int spacedim>
54 class TriangulationBase;
55}
56
57template <int dim, int spacedim>
59class DoFHandler;
60template <int dim, int spacedim, bool lda>
61class DoFCellAccessor;
62
63
64template <int dim, int spacedim>
65class Manifold;
66
67template <int dim, int spacedim>
68class Mapping;
69#endif
70
71namespace internal
72{
73 namespace TriangulationImplementation
74 {
75 class TriaObjects;
76 struct Implementation;
77 struct ImplementationMixedMesh;
78 } // namespace TriangulationImplementation
79
80 namespace TriaAccessorImplementation
81 {
82 struct Implementation;
83
89 template <int structdim, int dim>
91 {
92 struct type
93 {
97 type() = default;
98
102 type(const int level)
103 {
105 (void)level; // removes -Wunused-parameter warning in optimized mode
106 }
107
111 operator int() const
112 {
113 return 0;
114 }
115
116 void
118 {
120 }
121
122 void
124 {
126 }
127 };
128 };
129
130
136 template <int dim>
137 struct PresentLevelType<dim, dim>
138 {
139 using type = int;
140 };
141 } // namespace TriaAccessorImplementation
142} // namespace internal
143template <int structdim, int dim, int spacedim>
144class TriaAccessor;
145template <int dim, int spacedim>
146class TriaAccessor<0, dim, spacedim>;
147template <int spacedim>
148class TriaAccessor<0, 1, spacedim>;
149
150// note: the file tria_accessor.templates.h is included at the end of
151// this file. this includes a lot of templates. originally, this was
152// only done in debug mode, but led to cyclic reduction problems and
153// so is now on by default.
154
155
160{
165 "The operation you are attempting can only be performed for "
166 "(cell, face, or edge) iterators that point to valid "
167 "objects. These objects need not necessarily be active, "
168 "i.e., have no children, but they need to be part of a "
169 "triangulation. (The objects pointed to by an iterator "
170 "may -- after coarsening -- also be objects that used "
171 "to be part of a triangulation, but are now no longer "
172 "used. Their memory location may have been retained "
173 "for re-use upon the next mesh refinement, but is "
174 "currently unused.)");
185 "The operation you are attempting can only be performed for "
186 "(cell, face, or edge) iterators that point to 'active' "
187 "objects. 'Active' objects are those that do not have "
188 "children (in the case of cells), or that are part of "
189 "an active cell (in the case of faces or edges). However, "
190 "the object on which you are trying the current "
191 "operation is not 'active' in this sense.");
198 "The operation you are attempting can only be performed for "
199 "(cell, face, or edge) iterators that have children, "
200 "but the object on which you are trying the current "
201 "operation does not have any.");
209 "The operation you are attempting can only be performed for "
210 "(cell, face, or edge) iterators that have a parent object, "
211 "but the object on which you are trying the current "
212 "operation does not have one -- i.e., it is on the "
213 "coarsest level of the triangulation.");
218 int,
219 << "You can only set the child index if the cell does not "
220 << "currently have children registered; or you can clear it. "
221 << "The given index was " << arg1
222 << " (-1 means: clear children).");
226 template <typename AccessorType>
228 AccessorType,
229 << "You tried to dereference an iterator for which this "
230 << "is not possible. More information on this iterator: "
231 << "index=" << arg1.index() << ", state="
232 << (arg1.state() == IteratorState::valid ?
233 "valid" :
234 (arg1.state() == IteratorState::past_the_end ?
235 "past_the_end" :
236 "invalid")));
241 "Iterators can only be compared if they point to the same "
242 "triangulation, or if neither of them are associated "
243 "with a triangulation.");
244 // TODO: Write documentation!
249 // TODO: Write documentation!
269 // TODO: Write documentation!
275 int,
276 << "You can only set the child index of an even numbered child."
277 << "The number of the child given was " << arg1 << '.');
278} // namespace TriaAccessorExceptions
279
280
306template <int structdim, int dim, int spacedim = dim>
308{
309public:
316 static constexpr unsigned int space_dimension = spacedim;
317
323 static constexpr unsigned int dimension = dim;
324
330 static const unsigned int structure_dimension = structdim;
331
341 void
342 operator=(const TriaAccessorBase *) = delete;
343
344protected:
350 using AccessorData = void;
351
356 const int level = -1,
357 const int index = -1,
358 const AccessorData * = nullptr);
359
364
372 void
374
380
391 bool
392 operator<(const TriaAccessorBase &other) const;
393
394protected:
398 bool
400
404 bool
406
420 void
422
430 void
440 objects() const;
441
442public:
448 using LocalData = void *;
449
473 int
474 level() const;
475
502 int
503 index() const;
504
510 state() const;
511
518
522protected:
527 typename ::internal::TriaAccessorImplementation::
528 PresentLevelType<structdim, dim>::type present_level;
529
535
540
541private:
542 template <typename Accessor>
543 friend class TriaRawIterator;
544 template <typename Accessor>
545 friend class TriaIterator;
546 template <typename Accessor>
547 friend class TriaActiveIterator;
548};
549
550
551
572template <int structdim, int dim, int spacedim = dim>
574{
575public:
582 static constexpr unsigned int space_dimension = spacedim;
583
589 static constexpr unsigned int dimension = dim;
590
596 static const unsigned int structure_dimension = structdim;
597
603 using AccessorData = void;
604
612 InvalidAccessor(const void *parent = nullptr,
613 const int level = -1,
614 const int index = -1,
615 const AccessorData *local_data = nullptr);
616
625
630 template <typename OtherAccessor>
631 InvalidAccessor(const OtherAccessor &);
632
636 void
638
642 bool
644 bool
646
650 void
651 operator++() const;
652 void
653 operator--() const;
654
660 state();
661
662
667 static int
668 level();
669
674 static int
675 index();
676
681 bool
682 used() const;
683
688 bool
690
695 manifold_id() const;
696
700 unsigned int
701 user_index() const;
702
706 void
707 set_user_index(const unsigned int p) const;
708
712 void
714
719 vertex(const unsigned int i) const;
720
725 void *
726 line(const unsigned int i) const;
727
732 void *
733 quad(const unsigned int i) const;
734};
735
736
737
755template <int structdim, int dim, int spacedim>
756class TriaAccessor : public TriaAccessorBase<structdim, dim, spacedim>
757{
758public:
764
769 const int level = -1,
770 const int index = -1,
771 const AccessorData *local_data = nullptr);
772
777 TriaAccessor(const TriaAccessor &) = default;
778
782 TriaAccessor(TriaAccessor &&) = default; // NOLINT
783
796 template <int structdim2, int dim2, int spacedim2>
798
803 template <int structdim2, int dim2, int spacedim2>
805
816 operator=(const TriaAccessor &) = delete;
817
822 operator=(TriaAccessor &&) = default; // NOLINT
823
827 ~TriaAccessor() = default;
828
835 bool
836 used() const;
837
850 vertex_iterator(const unsigned int i) const;
851
867 unsigned int
868 vertex_index(const unsigned int i) const;
869
908 vertex(const unsigned int i) const;
909
913 typename ::internal::TriangulationImplementation::
914 Iterators<dim, spacedim>::line_iterator
915 line(const unsigned int i) const;
916
923 unsigned int
924 line_index(const unsigned int i) const;
925
929 typename ::internal::TriangulationImplementation::
930 Iterators<dim, spacedim>::quad_iterator
931 quad(const unsigned int i) const;
932
939 unsigned int
940 quad_index(const unsigned int i) const;
958 unsigned char
959 combined_face_orientation(const unsigned int face) const;
960
972 bool
973 face_orientation(const unsigned int face) const;
974
984 bool
985 face_flip(const unsigned int face) const;
986
996 bool
997 face_rotation(const unsigned int face) const;
998
1011 bool
1012 line_orientation(const unsigned int line) const;
1027 bool
1029
1034 unsigned int
1035 n_children() const;
1036
1050 unsigned int
1052
1066 unsigned int
1068
1073 child(const unsigned int i) const;
1074
1079 unsigned int
1082
1092 isotropic_child(const unsigned int i) const;
1093
1099
1105 int
1106 child_index(const unsigned int i) const;
1107
1113 int
1114 isotropic_child_index(const unsigned int i) const;
1138
1166 void
1168
1197 void
1199
1207 bool
1209
1221
1244
1262 void
1264
1278 void
1280
1297 bool
1299
1305 void
1307
1313 void
1315
1321 void
1323
1329 void
1331
1337 void
1339
1351 void
1352 set_user_pointer(void *p) const;
1353
1359 void
1361
1379 void *
1381
1403 void
1405
1412 void
1414
1424 void
1425 set_user_index(const unsigned int p) const;
1426
1432 void
1434
1446 unsigned int
1447 user_index() const;
1448
1466 void
1467 recursively_set_user_index(const unsigned int p) const;
1468
1477 void
1513 double
1514 diameter() const;
1515
1542 std::pair<Point<spacedim>, double>
1544
1555
1565 double
1566 extent_in_direction(const unsigned int axis) const;
1567
1571 double
1573
1588 intermediate_point(const Point<structdim> &coordinates) const;
1589
1614
1650 center(const bool respect_manifold = false,
1651 const bool interpolate_from_surrounding = false) const;
1652
1671 barycenter() const;
1672
1703 double
1704 measure() const;
1705
1720 bool
1723
1729
1733 unsigned int
1734 n_vertices() const;
1735
1739 unsigned int
1740 n_lines() const;
1741
1751 unsigned int
1752 n_faces() const;
1753
1760
1767
1776
1781private:
1786 void
1788
1796 void
1798 const std::initializer_list<int> &new_indices) const;
1799
1803 void
1805 const std::initializer_list<unsigned int> &new_indices) const;
1806
1814 void
1815 set_line_orientation(const unsigned int line, const bool orientation) const;
1816
1824 void
1825 set_combined_face_orientation(const unsigned int face,
1826 const unsigned char combined_orientation) const;
1827
1831 void
1833
1837 void
1839
1848 void
1850
1858 void
1860
1867 void
1868 set_children(const unsigned int i, const int index) const;
1869
1874 void
1876
1877private:
1878 friend class Triangulation<dim, spacedim>;
1879
1880 friend struct ::internal::TriangulationImplementation::Implementation;
1881 friend struct ::internal::TriangulationImplementation::
1882 ImplementationMixedMesh;
1883 friend struct ::internal::TriaAccessorImplementation::Implementation;
1884};
1885
1886
1887
1906template <int dim, int spacedim>
1907class TriaAccessor<0, dim, spacedim>
1908{
1909public:
1915 static constexpr unsigned int space_dimension = spacedim;
1916
1922 static constexpr unsigned int dimension = dim;
1923
1929 static const unsigned int structure_dimension = 0;
1930
1934 using AccessorData = void;
1935
1941 const unsigned int vertex_index);
1942
1949 const int level = 0,
1950 const int index = 0,
1951 const AccessorData * = nullptr);
1952
1956 template <int structdim2, int dim2, int spacedim2>
1958
1962 template <int structdim2, int dim2, int spacedim2>
1964
1969 state() const;
1970
1975 static int
1977
1982 int
1983 index() const;
1984
1991
2001 void
2003
2007 void
2012 bool
2013 operator==(const TriaAccessor &) const;
2014
2018 bool
2019 operator!=(const TriaAccessor &) const;
2020
2048 unsigned int
2049 vertex_index(const unsigned int i = 0) const;
2050
2057 vertex(const unsigned int i = 0) const;
2058
2063 typename ::internal::TriangulationImplementation::
2064 Iterators<dim, spacedim>::line_iterator static line(const unsigned int);
2065
2069 static unsigned int
2070 line_index(const unsigned int i);
2071
2075 static typename ::internal::TriangulationImplementation::
2076 Iterators<dim, spacedim>::quad_iterator
2077 quad(const unsigned int i);
2078
2082 static unsigned int
2083 quad_index(const unsigned int i);
2084
2100 double
2101 diameter() const;
2102
2110 double
2111 extent_in_direction(const unsigned int axis) const;
2112
2121 center(const bool respect_manifold = false,
2122 const bool interpolate_from_surrounding = false) const;
2123
2132 double
2133 measure() const;
2150 static unsigned char
2151 combined_face_orientation(const unsigned int face);
2152
2156 static bool
2157 face_orientation(const unsigned int face);
2158
2162 static bool
2163 face_flip(const unsigned int face);
2164
2168 static bool
2169 face_rotation(const unsigned int face);
2170
2174 static bool
2175 line_orientation(const unsigned int line);
2176
2191 static bool
2193
2198 static unsigned int
2200
2205 static unsigned int
2207
2211 static unsigned int
2213
2217 static unsigned int
2219
2224 child(const unsigned int);
2225
2230 isotropic_child(const unsigned int);
2231
2235 static RefinementCase<0>
2237
2241 static int
2242 child_index(const unsigned int i);
2243
2247 static int
2248 isotropic_child_index(const unsigned int i);
2256 bool
2257 used() const;
2258
2259protected:
2267 void
2269
2278 bool
2279 operator<(const TriaAccessor &other) const;
2280
2285
2290
2291private:
2292 template <typename Accessor>
2293 friend class TriaRawIterator;
2294 template <typename Accessor>
2295 friend class TriaIterator;
2296 template <typename Accessor>
2298};
2299
2300
2301
2318template <int spacedim>
2319class TriaAccessor<0, 1, spacedim>
2320{
2321public:
2327 static constexpr unsigned int space_dimension = spacedim;
2328
2334 static constexpr unsigned int dimension = 1;
2335
2341 static const unsigned int structure_dimension = 0;
2342
2346 using AccessorData = void;
2347
2353 {
2365 right_vertex
2367
2380 const VertexKind vertex_kind,
2381 const unsigned int vertex_index);
2382
2389 const int = 0,
2390 const int = 0,
2391 const AccessorData * = nullptr);
2392
2396 template <int structdim2, int dim2, int spacedim2>
2398
2402 template <int structdim2, int dim2, int spacedim2>
2404
2409 void
2411
2417 void
2419
2427
2432 static int
2434
2439 int
2440 index() const;
2441
2448
2459 void
2460 operator++() const;
2461
2466 void
2467 operator--() const;
2471 bool
2472 operator==(const TriaAccessor &) const;
2473
2477 bool
2478 operator!=(const TriaAccessor &) const;
2479
2488 bool
2489 operator<(const TriaAccessor &other) const;
2490
2517 unsigned int
2518 vertex_index(const unsigned int i = 0) const;
2519
2526 vertex(const unsigned int i = 0) const;
2527
2533 center() const;
2534
2539 typename ::internal::TriangulationImplementation::
2540 Iterators<1, spacedim>::line_iterator static line(const unsigned int);
2541
2548 static unsigned int
2549 line_index(const unsigned int i);
2550
2554 static typename ::internal::TriangulationImplementation::
2555 Iterators<1, spacedim>::quad_iterator
2556 quad(const unsigned int i);
2557
2564 static unsigned int
2565 quad_index(const unsigned int i);
2566
2576 bool
2578
2595
2599 const Manifold<1, spacedim> &
2601
2610
2611
2623 bool
2624 user_flag_set() const;
2625
2631 void
2632 set_user_flag() const;
2633
2639 void
2640 clear_user_flag() const;
2641
2647 void
2649
2655 void
2657
2663 void
2664 clear_user_data() const;
2665
2677 void
2678 set_user_pointer(void *p) const;
2679
2685 void
2686 clear_user_pointer() const;
2687
2703 void *
2704 user_pointer() const;
2705
2727 void
2728 recursively_set_user_pointer(void *p) const;
2729
2736 void
2738
2748 void
2749 set_user_index(const unsigned int p) const;
2750
2756 void
2757 clear_user_index() const;
2758
2770 unsigned int
2771 user_index() const;
2772
2790 void
2791 recursively_set_user_index(const unsigned int p) const;
2792
2801 void
2819 static unsigned char
2820 combined_face_orientation(const unsigned int face);
2821
2825 static bool
2826 face_orientation(const unsigned int face);
2827
2831 static bool
2832 face_flip(const unsigned int face);
2833
2837 static bool
2838 face_rotation(const unsigned int face);
2839
2843 static bool
2844 line_orientation(const unsigned int line);
2845
2860 static bool
2862
2867 static unsigned int
2869
2874 static unsigned int
2876
2880 static unsigned int
2882
2886 static unsigned int
2888
2893 child(const unsigned int);
2894
2899 isotropic_child(const unsigned int);
2900
2904 static RefinementCase<0>
2906
2910 static int
2911 child_index(const unsigned int i);
2912
2916 static int
2917 isotropic_child_index(const unsigned int i);
2948 void
2950
2957 void
2959
2969 void
2971
2983 void
2992 bool
2993 used() const;
2994
3000
3004 unsigned int
3005 n_vertices() const;
3006
3010 unsigned int
3011 n_lines() const;
3012
3019
3026
3027protected:
3032
3038
3043};
3044
3045
3046
3062template <int dim, int spacedim = dim>
3063class CellAccessor : public TriaAccessor<dim, dim, spacedim>
3064{
3065public:
3070
3075
3087 const int level = -1,
3088 const int index = -1,
3089 const AccessorData *local_data = nullptr);
3090
3095
3108 template <int structdim2, int dim2, int spacedim2>
3110
3115 template <int structdim2, int dim2, int spacedim2>
3117
3122
3126 // NOLINTNEXTLINE OSX does not compile with noexcept
3128
3132 ~CellAccessor() = default;
3133
3145
3149 // NOLINTNEXTLINE OSX does not compile with noexcept
3152
3176 as_dof_handler_iterator(const DoFHandler<dim, spacedim> &dof_handler) const;
3177
3188 const DoFHandler<dim, spacedim> &dof_handler) const;
3189
3190
3207 child(const unsigned int i) const;
3208
3212 boost::container::small_vector<TriaIterator<CellAccessor<dim, spacedim>>,
3215
3219 TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>
3220 face(const unsigned int i) const;
3221
3226 unsigned int
3229
3233 boost::container::small_vector<
3234 TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>,
3237
3247 unsigned int
3248 face_index(const unsigned int i) const;
3249
3298 neighbor_child_on_subface(const unsigned int face_no,
3299 const unsigned int subface_no) const;
3300
3349 neighbor(const unsigned int face_no) const;
3350
3358 int
3359 neighbor_index(const unsigned int face_no) const;
3360
3368 int
3369 neighbor_level(const unsigned int face_no) const;
3370
3382 unsigned int
3383 neighbor_of_neighbor(const unsigned int face_no) const;
3384
3395 bool
3396 neighbor_is_coarser(const unsigned int face_no) const;
3397
3412 std::pair<unsigned int, unsigned int>
3413 neighbor_of_coarser_neighbor(const unsigned int neighbor) const;
3414
3421 unsigned int
3422 neighbor_face_no(const unsigned int neighbor) const;
3423
3427 static bool
3429
3443 bool
3444 has_periodic_neighbor(const unsigned int i) const;
3445
3463 periodic_neighbor(const unsigned int i) const;
3464
3473 neighbor_or_periodic_neighbor(const unsigned int i) const;
3474
3490 periodic_neighbor_child_on_subface(const unsigned int face_no,
3491 const unsigned int subface_no) const;
3492
3503 std::pair<unsigned int, unsigned int>
3504 periodic_neighbor_of_coarser_periodic_neighbor(const unsigned face_no) const;
3505
3511 int
3512 periodic_neighbor_index(const unsigned int i) const;
3513
3519 int
3520 periodic_neighbor_level(const unsigned int i) const;
3521
3536 unsigned int
3537 periodic_neighbor_of_periodic_neighbor(const unsigned int i) const;
3538
3544 unsigned int
3545 periodic_neighbor_face_no(const unsigned int i) const;
3546
3553 bool
3554 periodic_neighbor_is_coarser(const unsigned int i) const;
3555
3572 bool
3573 at_boundary(const unsigned int i) const;
3574
3583 bool
3584 at_boundary() const;
3585
3593 bool
3594 has_boundary_lines() const;
3622
3640 void
3643
3647 void
3649
3654 std::uint8_t
3656
3661 void
3662 set_refine_choice(const std::uint8_t refinement_choice = static_cast<char>(
3664
3668 void
3670
3678 bool
3680 const unsigned int face_no,
3681 const RefinementCase<dim - 1> &face_refinement_case =
3683
3689 bool
3690 flag_for_line_refinement(const unsigned int line_no) const;
3691
3701 subface_case(const unsigned int face_no) const;
3702
3706 bool
3708
3713 void
3715
3719 void
3744 material_id() const;
3745
3757 void
3758 set_material_id(const types::material_id new_material_id) const;
3759
3768 void
3769 recursively_set_material_id(const types::material_id new_material_id) const;
3797
3813 void
3814 set_subdomain_id(const types::subdomain_id new_subdomain_id) const;
3815
3826
3831 void
3833 const types::subdomain_id new_level_subdomain_id) const;
3834
3835
3851 void
3853 const types::subdomain_id new_subdomain_id) const;
3877
3887
3904 bool
3905 direction_flag() const;
3906
3932 unsigned int
3934
3942 int
3943 parent_index() const;
3944
3951 parent() const;
3952
3972 bool
3973 is_active() const;
3974
3994 bool
3996
4001 bool
4003
4037 bool
4038 is_ghost() const;
4039
4045 bool
4047
4074 bool
4076
4083 bool
4085
4099 bool
4101
4110 void
4111 set_neighbor(const unsigned int i,
4112 const TriaIterator<CellAccessor<dim, spacedim>> &pointer) const;
4113
4127 CellId
4128 id() const;
4129
4130 using TriaAccessor<dim, dim, spacedim>::diameter;
4131
4135 double
4136 diameter(const Mapping<dim, spacedim> &mapping) const;
4137
4155
4156protected:
4172 unsigned int
4173 neighbor_of_neighbor_internal(const unsigned int neighbor) const;
4174
4180 template <int dim_, int spacedim_>
4181 bool
4182 point_inside_codim(const Point<spacedim_> &p) const;
4183
4184
4185
4186private:
4191 void
4192 set_active_cell_index(const unsigned int active_cell_index) const;
4193
4197 void
4199
4203 void
4205
4209 void
4210 set_parent(const unsigned int parent_index);
4211
4222 void
4223 set_direction_flag(const bool new_direction_flag) const;
4224
4225 friend class Triangulation<dim, spacedim>;
4226
4227 friend class parallel::TriangulationBase<dim, spacedim>;
4228
4229 friend struct ::internal::TriangulationImplementation::Implementation;
4230 friend struct ::internal::TriangulationImplementation::
4231 ImplementationMixedMesh;
4232};
4233
4234
4235
4236/* ----- declaration of explicit specializations and general templates ----- */
4237
4238
4239template <int structdim, int dim, int spacedim>
4240template <typename OtherAccessor>
4242 const OtherAccessor &)
4243{
4244 Assert(false,
4245 ExcMessage("You are attempting an illegal conversion between "
4246 "iterator/accessor types. The constructor you call "
4247 "only exists to make certain template constructs "
4248 "easier to write as dimension independent code but "
4249 "the conversion is not valid in the current context."));
4250}
4251
4252
4253
4254template <int structdim, int dim, int spacedim>
4255template <int structdim2, int dim2, int spacedim2>
4258{
4259 Assert(false,
4260 ExcMessage("You are attempting an illegal conversion between "
4261 "iterator/accessor types. The constructor you call "
4262 "only exists to make certain template constructs "
4263 "easier to write as dimension independent code but "
4264 "the conversion is not valid in the current context."));
4265}
4266
4267
4268
4269template <int dim, int spacedim>
4270template <int structdim2, int dim2, int spacedim2>
4273{
4274 Assert(false,
4275 ExcMessage("You are attempting an illegal conversion between "
4276 "iterator/accessor types. The constructor you call "
4277 "only exists to make certain template constructs "
4278 "easier to write as dimension independent code but "
4279 "the conversion is not valid in the current context."));
4280}
4281
4282
4283
4284template <int structdim, int dim, int spacedim>
4285template <int structdim2, int dim2, int spacedim2>
4288{
4289 Assert(false,
4290 ExcMessage("You are attempting an illegal conversion between "
4291 "iterator/accessor types. The constructor you call "
4292 "only exists to make certain template constructs "
4293 "easier to write as dimension independent code but "
4294 "the conversion is not valid in the current context."));
4295}
4296
4297
4298
4299template <int dim, int spacedim>
4300template <int structdim2, int dim2, int spacedim2>
4303{
4304 Assert(false,
4305 ExcMessage("You are attempting an illegal conversion between "
4306 "iterator/accessor types. The constructor you call "
4307 "only exists to make certain template constructs "
4308 "easier to write as dimension independent code but "
4309 "the conversion is not valid in the current context."));
4310}
4311
4312
4313#ifndef DOXYGEN
4314
4315template <>
4316bool
4318template <>
4319bool
4321template <>
4322bool
4324template <>
4325bool
4327template <>
4328bool
4330template <>
4331bool
4333// -------------------------------------------------------------------
4334
4335template <>
4336void
4338
4339#endif // DOXYGEN
4340
4342
4343// include more templates in debug and optimized mode
4344#include <deal.II/grid/tria_accessor.templates.h>
4345
4346#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
std::uint8_t refine_choice() 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
TriaActiveIterator< DoFCellAccessor< dim, spacedim, false > > as_dof_handler_iterator(const DoFHandler< dim, spacedim > &dof_handler) 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
void clear_refine_choice() 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
void set_refine_choice(const std::uint8_t refinement_choice=static_cast< char >(IsotropicRefinementChoice::isotropic_refinement)) 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
TriaIterator< DoFCellAccessor< dim, spacedim, true > > as_dof_handler_level_iterator(const DoFHandler< dim, spacedim > &dof_handler) 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
void set_manifold_id(const types::manifold_id) const
InvalidAccessor(const InvalidAccessor &)
static constexpr unsigned int space_dimension
static int level()
bool operator==(const InvalidAccessor &) const
Point< spacedim > & vertex(const unsigned int i) const
void * quad(const unsigned int i) const
void operator++() const
void * line(const unsigned int i) const
static IteratorState::IteratorStates state()
bool operator!=(const InvalidAccessor &) const
bool has_children() const
unsigned int user_index() const
static const unsigned int structure_dimension
void copy_from(const InvalidAccessor &)
void set_user_index(const unsigned int p) const
static int index()
types::manifold_id manifold_id() const
InvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
static constexpr unsigned int dimension
void operator--() const
bool used() const
InvalidAccessor(const OtherAccessor &)
Abstract base class for mapping classes.
Definition mapping.h:318
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 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 > &)
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 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()
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 > &)
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 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 unsigned int max_refinement_depth()
static unsigned int quad_index(const unsigned int i)
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
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
void clear_children() const
void set_boundary_id_internal(const types::boundary_id id) const
void set_user_index(const unsigned int p) const
TriaAccessor(TriaAccessor &&)=default
void clear_user_pointer() const
unsigned int n_active_descendants() 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
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 set_bounding_object_indices(const std::initializer_list< unsigned int > &new_indices) 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
Point< spacedim > center(const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
~TriaAccessor()=default
void set_refinement_case(const RefinementCase< structdim > &ref_case) const
void clear_user_index() const
double minimum_vertex_distance() const
double measure() const
void set_all_boundary_ids(const types::boundary_id) const
void set_bounding_object_indices(const std::initializer_list< int > &new_indices) 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
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_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:175
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
unsigned int level
Definition grid_out.cc:4626
static ::ExceptionBase & ExcCantCompareIterators()
#define DeclException0(Exception0)
Definition exceptions.h:466
static ::ExceptionBase & ExcCellHasNoParent()
static ::ExceptionBase & ExcCellFlaggedForRefinement()
static ::ExceptionBase & ExcCellHasNoChildren()
static ::ExceptionBase & ExcCellNotUsed()
static ::ExceptionBase & ExcNoPeriodicNeighbor()
static ::ExceptionBase & ExcCantSetChildren(int arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcNeighborIsNotCoarser()
static ::ExceptionBase & ExcSetOnlyEvenChildren(int arg1)
static ::ExceptionBase & ExcNeighborIsCoarser()
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:489
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcFacesHaveNoLevel()
static ::ExceptionBase & ExcDereferenceInvalidObject(AccessorType arg1)
static ::ExceptionBase & ExcCellFlaggedForCoarsening()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:511
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
bool line_orientation(const unsigned int line) const
static bool face_orientation(const unsigned int face)
Always return false.
unsigned char combined_face_orientation(const unsigned int face) const
bool face_rotation(const unsigned int face) const
bool face_orientation(const unsigned int face) const
static bool face_flip(const unsigned int face)
Always return false.
static bool line_orientation(const unsigned int line)
Always return false.
static unsigned char combined_face_orientation(const unsigned int face)
Always return 0.
static bool face_flip(const unsigned int face)
Always return false.
static bool line_orientation(const unsigned int line)
Always return false.
static bool face_rotation(const unsigned int face)
Always return false.
static bool face_rotation(const unsigned int face)
Always return false.
static bool face_orientation(const unsigned int face)
Always return false.
static unsigned char combined_face_orientation(const unsigned int face)
Always return 0.
bool face_flip(const unsigned int face) const
void set_combined_face_orientation(const unsigned int face, const unsigned char combined_orientation) const
#define DEAL_II_ASSERT_UNREACHABLE()
@ valid
Iterator points to a valid object.
boost::integer_range< IncrementableType > iota_view
Definition iota_view.h:45