deal.II version GIT relicensing-2287-g6548a49e0a 2024-12-20 18:30: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.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_h
16#define dealii_tria_h
17
18
19#include <deal.II/base/config.h>
20
26#include <deal.II/base/point.h>
27
33
34#include <boost/serialization/map.hpp>
35#include <boost/serialization/split_member.hpp>
36#include <boost/serialization/unique_ptr.hpp>
37#include <boost/serialization/vector.hpp>
38#include <boost/signals2.hpp>
39
40#include <bitset>
41#include <functional>
42#include <list>
43#include <map>
44#include <memory>
45#include <numeric>
46#include <vector>
47
48
50
51#ifdef signals
52# error \
53 "The name 'signals' is already defined. You are most likely using the QT library \
54and using the 'signals' keyword. You can either #include the Qt headers (or any conflicting headers) \
55*after* the deal.II headers or you can define the 'QT_NO_KEYWORDS' macro and use the 'Q_SIGNALS' macro."
56#endif
57
58// Forward declarations
59#ifndef DOXYGEN
60template <int dim, int spacedim>
61class Manifold;
62
63template <int dim>
64struct CellData;
65
66struct SubCellData;
67
69{
70 template <int, int>
71 struct Description;
72}
73
74namespace GridTools
75{
76 template <typename CellIterator>
77 struct PeriodicFacePair;
78}
79
80template <int, int, int>
81class TriaAccessor;
82template <int spacedim>
83class TriaAccessor<0, 1, spacedim>;
84template <int, int, int>
86
87namespace internal
88{
89 namespace TriangulationImplementation
90 {
91 class TriaFaces;
92
93 class TriaObjects;
94
95 template <int, int>
96 class Policy;
97
103 struct Implementation;
104 struct ImplementationMixedMesh;
105 } // namespace TriangulationImplementation
106
107 namespace TriaAccessorImplementation
108 {
109 struct Implementation;
110 }
111} // namespace internal
112#endif
113
114
115/*------------------------------------------------------------------------*/
116
117
118namespace internal
119{
124 namespace TriangulationImplementation
125 {
137 template <int dim>
139 {};
140
152 template <>
153 struct NumberCache<1>
154 {
158 unsigned int n_levels;
159
163 unsigned int n_lines;
164
168 std::vector<unsigned int> n_lines_level;
169
173 unsigned int n_active_lines;
174
178 std::vector<unsigned int> n_active_lines_level;
179
183 std::shared_ptr<const Utilities::MPI::Partitioner>
185
189 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
191
195 NumberCache();
196
201 std::size_t
202 memory_consumption() const;
203
209 template <class Archive>
210 void
211 serialize(Archive &ar, const unsigned int version);
212 };
213
214
227 template <>
228 struct NumberCache<2> : public NumberCache<1>
229 {
233 unsigned int n_quads;
234
238 std::vector<unsigned int> n_quads_level;
239
243 unsigned int n_active_quads;
244
248 std::vector<unsigned int> n_active_quads_level;
249
253 NumberCache();
254
259 std::size_t
260 memory_consumption() const;
261
267 template <class Archive>
268 void
269 serialize(Archive &ar, const unsigned int version);
270 };
271
272
286 template <>
287 struct NumberCache<3> : public NumberCache<2>
288 {
292 unsigned int n_hexes;
293
297 std::vector<unsigned int> n_hexes_level;
298
302 unsigned int n_active_hexes;
303
307 std::vector<unsigned int> n_active_hexes_level;
308
312 NumberCache();
313
318 std::size_t
319 memory_consumption() const;
320
326 template <class Archive>
327 void
328 serialize(Archive &ar, const unsigned int version);
329 };
330 } // namespace TriangulationImplementation
331
332
336 template <int dim, int spacedim = dim>
339 {
341
347
353
355 std::function<std::vector<char>(cell_iterator, CellStatus)>;
356
361 std::vector<pack_callback_t> pack_callbacks_fixed;
362 std::vector<pack_callback_t> pack_callbacks_variable;
363 };
364
375 template <int dim, int spacedim = dim>
378 {
379 public:
381
386 static inline constexpr unsigned int version_number = 5;
387
393 using cell_relation_t = typename std::pair<cell_iterator, CellStatus>;
394
396
405 void
406 pack_data(
407 const std::vector<cell_relation_t> &cell_relations,
408 const std::vector<
410 &pack_callbacks_fixed,
411 const std::vector<
413 &pack_callbacks_variable,
414 const MPI_Comm &mpi_communicator);
415
423 void
424 unpack_cell_status(std::vector<cell_relation_t> &cell_relations) const;
425
438 void
439 unpack_data(
440 const std::vector<cell_relation_t> &cell_relations,
441 const unsigned int handle,
442 const std::function<
443 void(const cell_iterator &,
444 const CellStatus &,
445 const boost::iterator_range<std::vector<char>::const_iterator> &)>
446 &unpack_callback) const;
447
462 void
463 save(const unsigned int global_first_cell,
464 const unsigned int global_num_cells,
465 const std::string &file_basename,
466 const MPI_Comm &mpi_communicator) const;
467
486 void
487 load(const unsigned int global_first_cell,
488 const unsigned int global_num_cells,
489 const unsigned int local_num_cells,
490 const std::string &file_basename,
491 const unsigned int n_attached_deserialize_fixed,
492 const unsigned int n_attached_deserialize_variable,
493 const MPI_Comm &mpi_communicator);
494
501 void
502 clear();
503
508
519 std::vector<unsigned int> sizes_fixed_cumulative;
520
525 std::vector<char> src_data_fixed;
526 std::vector<char> dest_data_fixed;
527
532 std::vector<int> src_sizes_variable;
533 std::vector<int> dest_sizes_variable;
534 std::vector<char> src_data_variable;
535 std::vector<char> dest_data_variable;
536 };
537} // namespace internal
538
539
540/*------------------------------------------------------------------------*/
541
542
1320template <int dim, int spacedim = dim>
1323{
1324private:
1331
1332public:
1338 {
1343 none = 0x0,
1387 limit_level_difference_at_vertices = 0x1,
1408 eliminate_unrefined_islands = 0x2,
1424 patch_level_1 = 0x4,
1448 coarsest_level_1 = 0x8,
1473 allow_anisotropic_smoothing = 0x10,
1506 eliminate_refined_inner_islands = 0x100,
1511 eliminate_refined_boundary_islands = 0x200,
1517 do_not_produce_unrefined_islands = 0x400,
1518
1523 smoothing_on_refinement =
1524 (limit_level_difference_at_vertices | eliminate_unrefined_islands),
1529 smoothing_on_coarsening =
1530 (eliminate_refined_inner_islands | eliminate_refined_boundary_islands |
1531 do_not_produce_unrefined_islands),
1532
1538 maximum_smoothing = 0xffff ^ allow_anisotropic_smoothing
1540
1557
1563
1581
1595 using face_iterator = TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1596
1609 TriaActiveIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1610
1620
1635
1643 using line_iterator = typename IteratorSelector::line_iterator;
1644
1658 using active_line_iterator = typename IteratorSelector::active_line_iterator;
1659
1667 using quad_iterator = typename IteratorSelector::quad_iterator;
1668
1682 using active_quad_iterator = typename IteratorSelector::active_quad_iterator;
1683
1691 using hex_iterator = typename IteratorSelector::hex_iterator;
1692
1702 using active_hex_iterator = typename IteratorSelector::active_hex_iterator;
1703
1723 {
1730 virtual ~DistortedCellList() noexcept override;
1731
1736 std::list<typename Triangulation<dim, spacedim>::cell_iterator>
1737 distorted_cells;
1738 };
1739
1743 static constexpr unsigned int dimension = dim;
1744
1748 static constexpr unsigned int space_dimension = spacedim;
1749
1764 Triangulation(const MeshSmoothing smooth_grid = none,
1765 const bool check_for_distorted_cells = false);
1766
1782 Triangulation(const Triangulation<dim, spacedim> &) = delete;
1783
1790 Triangulation(Triangulation<dim, spacedim> &&tria) noexcept;
1791
1796 operator=(Triangulation<dim, spacedim> &&tria) noexcept;
1797
1801 virtual ~Triangulation() override;
1802
1809 virtual void
1810 clear();
1811
1816 virtual MPI_Comm
1817 get_mpi_communicator() const;
1818
1825 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
1826 "Access the MPI communicator with get_mpi_communicator() instead.")
1827 MPI_Comm
1828 get_communicator() const;
1829
1835 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1836 global_active_cell_index_partitioner() const;
1837
1843 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1844 global_level_cell_index_partitioner(const unsigned int level) const;
1845
1850 virtual void
1851 set_mesh_smoothing(const MeshSmoothing mesh_smoothing);
1852
1856 virtual const MeshSmoothing &
1857 get_mesh_smoothing() const;
1858
1881 void
1882 set_manifold(const types::manifold_id number,
1883 const Manifold<dim, spacedim> &manifold_object);
1884
1901 void
1902 reset_manifold(const types::manifold_id manifold_number);
1903
1920 void
1921 reset_all_manifolds();
1922
1931 void
1932 set_all_manifold_ids(const types::manifold_id number);
1933
1942 void
1943 set_all_manifold_ids_on_boundary(const types::manifold_id number);
1944
1954 void
1955 set_all_manifold_ids_on_boundary(const types::boundary_id b_id,
1956 const types::manifold_id number);
1957
1971 const Manifold<dim, spacedim> &
1972 get_manifold(const types::manifold_id number) const;
1973
1984 virtual std::vector<types::boundary_id>
1985 get_boundary_ids() const;
1986
1999 virtual std::vector<types::manifold_id>
2000 get_manifold_ids() const;
2001
2023 virtual void
2024 copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
2025
2074 virtual void
2075 create_triangulation(const std::vector<Point<spacedim>> &vertices,
2076 const std::vector<CellData<dim>> &cells,
2077 const SubCellData &subcelldata);
2078
2091 virtual void
2092 create_triangulation(
2093 const TriangulationDescription::Description<dim, spacedim>
2094 &construction_data);
2095
2104 void
2105 flip_all_direction_flags();
2106
2118 void
2119 set_all_refine_flags();
2120
2143 void
2144 refine_global(const unsigned int times = 1);
2145
2158 void
2159 coarsen_global(const unsigned int times = 1);
2160
2192 virtual void
2193 execute_coarsening_and_refinement();
2194
2224 virtual bool
2225 prepare_coarsening_and_refinement();
2226
2242
2247 static constexpr auto CELL_PERSIST DEAL_II_DEPRECATED =
2249
2254 static constexpr auto CELL_REFINE DEAL_II_DEPRECATED =
2256
2261 static constexpr auto CELL_COARSEN DEAL_II_DEPRECATED =
2263
2268 static constexpr auto CELL_INVALID DEAL_II_DEPRECATED =
2270
2271
2277 template <typename T>
2279 {
2280 using result_type = T;
2281
2282 template <typename InputIterator>
2283 T
2284 operator()(InputIterator first, InputIterator last) const
2285 {
2286 return std::accumulate(first, last, T());
2287 }
2288 };
2289
2299 struct Signals
2300 {
2308 boost::signals2::signal<void()> create;
2309
2317 boost::signals2::signal<void()> pre_refinement;
2318
2324 boost::signals2::signal<void()> post_refinement;
2325
2332 boost::signals2::signal<void()> pre_partition;
2333
2341 boost::signals2::signal<void()> mesh_movement;
2342
2350 boost::signals2::signal<void(
2351 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2353
2360 boost::signals2::signal<void(
2361 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2363
2370 boost::signals2::signal<void(
2371 const Triangulation<dim, spacedim> &destination_tria)>
2373
2387 boost::signals2::signal<void()> clear;
2388
2398 boost::signals2::signal<void()> any_change;
2399
2436 boost::signals2::signal<unsigned int(const cell_iterator &,
2437 const ::CellStatus),
2440
2452 boost::signals2::signal<void()> pre_distributed_refinement;
2453
2461 boost::signals2::signal<void()> post_p4est_refinement;
2462
2471 boost::signals2::signal<void()> post_distributed_refinement;
2472
2478 boost::signals2::signal<void()> pre_distributed_repartition;
2479
2485 boost::signals2::signal<void()> post_distributed_repartition;
2486
2493 boost::signals2::signal<void()> pre_distributed_save;
2494
2500 boost::signals2::signal<void()> post_distributed_save;
2501
2508 boost::signals2::signal<void()> pre_distributed_load;
2509
2515 boost::signals2::signal<void()> post_distributed_load;
2516 };
2517
2527
2539 void
2540 save_refine_flags(std::ostream &out) const;
2541
2545 void
2546 save_refine_flags(std::vector<bool> &v) const;
2547
2551 void
2552 load_refine_flags(std::istream &in);
2553
2557 void
2558 load_refine_flags(const std::vector<bool> &v);
2559
2563 void
2564 save_coarsen_flags(std::ostream &out) const;
2565
2569 void
2570 save_coarsen_flags(std::vector<bool> &v) const;
2571
2575 void
2576 load_coarsen_flags(std::istream &out);
2577
2581 void
2582 load_coarsen_flags(const std::vector<bool> &v);
2583
2588 bool
2590
2602 void
2604
2610 void
2611 save_user_flags(std::ostream &out) const;
2612
2618 void
2619 save_user_flags(std::vector<bool> &v) const;
2620
2625 void
2626 load_user_flags(std::istream &in);
2627
2632 void
2633 load_user_flags(const std::vector<bool> &v);
2634
2639 void
2641
2646 void
2647 save_user_flags_line(std::ostream &out) const;
2648
2654 void
2655 save_user_flags_line(std::vector<bool> &v) const;
2656
2661 void
2662 load_user_flags_line(std::istream &in);
2663
2668 void
2669 load_user_flags_line(const std::vector<bool> &v);
2670
2675 void
2677
2682 void
2683 save_user_flags_quad(std::ostream &out) const;
2684
2690 void
2691 save_user_flags_quad(std::vector<bool> &v) const;
2692
2697 void
2698 load_user_flags_quad(std::istream &in);
2699
2704 void
2705 load_user_flags_quad(const std::vector<bool> &v);
2706
2707
2712 void
2714
2719 void
2720 save_user_flags_hex(std::ostream &out) const;
2721
2727 void
2728 save_user_flags_hex(std::vector<bool> &v) const;
2729
2734 void
2735 load_user_flags_hex(std::istream &in);
2736
2741 void
2742 load_user_flags_hex(const std::vector<bool> &v);
2743
2749 void
2751
2757 void
2758 save_user_indices(std::vector<unsigned int> &v) const;
2759
2764 void
2765 load_user_indices(const std::vector<unsigned int> &v);
2766
2772 void
2773 save_user_pointers(std::vector<void *> &v) const;
2774
2779 void
2780 load_user_pointers(const std::vector<void *> &v);
2781
2787 void
2788 save_user_indices_line(std::vector<unsigned int> &v) const;
2789
2794 void
2795 load_user_indices_line(const std::vector<unsigned int> &v);
2796
2802 void
2803 save_user_indices_quad(std::vector<unsigned int> &v) const;
2804
2809 void
2810 load_user_indices_quad(const std::vector<unsigned int> &v);
2811
2817 void
2818 save_user_indices_hex(std::vector<unsigned int> &v) const;
2819
2824 void
2825 load_user_indices_hex(const std::vector<unsigned int> &v);
2831 void
2832 save_user_pointers_line(std::vector<void *> &v) const;
2833
2838 void
2839 load_user_pointers_line(const std::vector<void *> &v);
2840
2846 void
2847 save_user_pointers_quad(std::vector<void *> &v) const;
2848
2853 void
2854 load_user_pointers_quad(const std::vector<void *> &v);
2855
2861 void
2862 save_user_pointers_hex(std::vector<void *> &v) const;
2863
2868 void
2869 load_user_pointers_hex(const std::vector<void *> &v);
2870
2896 begin(const unsigned int level = 0) const;
2897
2929 begin_active(const unsigned int level = 0) const;
2930
2936 end() const;
2937
2957 end(const unsigned int level) const;
2958
2979 end_active(const unsigned int level) const;
2980
2981
2986 last() const;
2987
2993
3002 create_cell_iterator(const CellId &cell_id) const;
3003
3014 bool
3015 contains_cell(const CellId &cell_id) const;
3035
3073
3090 cell_iterators_on_level(const unsigned int level) const;
3091
3108 active_cell_iterators_on_level(const unsigned int level) const;
3109
3112 /*-------------------------------------------------------------------------*/
3113
3123 begin_face() const;
3124
3130
3136 end_face() const;
3137
3158
3161 /*-------------------------------------------------------------------------*/
3162
3174
3182
3189 end_vertex() const;
3190
3210 unsigned int
3211 n_lines() const;
3212
3216 unsigned int
3217 n_lines(const unsigned int level) const;
3218
3222 unsigned int
3224
3228 unsigned int
3229 n_active_lines(const unsigned int level) const;
3230
3234 unsigned int
3235 n_quads() const;
3236
3240 unsigned int
3241 n_quads(const unsigned int level) const;
3242
3246 unsigned int
3248
3252 unsigned int
3253 n_active_quads(const unsigned int level) const;
3254
3258 unsigned int
3259 n_hexs() const;
3260
3265 unsigned int
3266 n_hexs(const unsigned int level) const;
3267
3271 unsigned int
3273
3278 unsigned int
3279 n_active_hexs(const unsigned int level) const;
3280
3285 unsigned int
3286 n_cells() const;
3287
3292 unsigned int
3293 n_cells(const unsigned int level) const;
3294
3299 unsigned int
3301
3311
3312
3317 unsigned int
3318 n_active_cells(const unsigned int level) const;
3319
3324 virtual types::coarse_cell_id
3326
3332 unsigned int
3333 n_faces() const;
3334
3340 unsigned int
3342
3360 unsigned int
3361 n_levels() const;
3362
3369 virtual unsigned int
3371
3381 virtual bool
3383
3391 unsigned int
3392 n_vertices() const;
3393
3402 const std::vector<Point<spacedim>> &
3404
3409 unsigned int
3411
3415 bool
3416 vertex_used(const unsigned int index) const;
3417
3422 const std::vector<bool> &
3424
3436 unsigned int
3438
3445 virtual types::subdomain_id
3447
3459
3466
3467
3484 unsigned int
3486
3496 unsigned int
3497 n_raw_lines(const unsigned int level) const;
3498
3508 unsigned int
3510
3520 unsigned int
3521 n_raw_quads(const unsigned int level) const;
3522
3532 unsigned int
3533 n_raw_hexs(const unsigned int level) const;
3534
3544 unsigned int
3545 n_raw_cells(const unsigned int level) const;
3546
3558 unsigned int
3560
3571 virtual std::size_t
3573
3583 template <class Archive>
3584 void
3585 save(Archive &ar, const unsigned int version) const;
3586
3604 template <class Archive>
3605 void
3606 load(Archive &ar, const unsigned int version);
3607
3608
3617 virtual void
3618 save(const std::string &file_basename) const;
3619
3623 virtual void
3624 load(const std::string &file_basename);
3625
3626
3642 virtual void
3644 const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3645
3649 const std::map<
3650 std::pair<cell_iterator, unsigned int>,
3651 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>> &
3653
3658 const std::vector<ReferenceCell> &
3660
3665 bool
3667
3672 bool
3674
3680 bool
3682
3683#ifdef DOXYGEN
3691 template <class Archive>
3692 void
3693 serialize(Archive &archive, const unsigned int version);
3694#else
3695 // This macro defines the serialize() method that is compatible with
3696 // the templated save() and load() method that have been implemented.
3697 BOOST_SERIALIZATION_SPLIT_MEMBER()
3698#endif
3699
3704public:
3812 unsigned int
3814 const std::function<std::vector<char>(const cell_iterator &,
3815 const ::CellStatus)>
3816 &pack_callback,
3817 const bool returns_variable_size_data);
3818
3868 void
3870 const unsigned int handle,
3871 const std::function<
3872 void(const cell_iterator &,
3873 const ::CellStatus,
3874 const boost::iterator_range<std::vector<char>::const_iterator> &)>
3875 &unpack_callback);
3876
3878
3879protected:
3887 void
3888 save_attached_data(const unsigned int global_first_cell,
3889 const unsigned int global_num_cells,
3890 const std::string &file_basename) const;
3891
3900 void
3901 load_attached_data(const unsigned int global_first_cell,
3902 const unsigned int global_num_cells,
3903 const unsigned int local_num_cells,
3904 const std::string &file_basename,
3905 const unsigned int n_attached_deserialize_fixed,
3906 const unsigned int n_attached_deserialize_variable);
3907
3917 void
3919
3925 void
3927
3928
3933 void
3935
3941 std::vector<typename internal::CellAttachedDataSerializer<dim, spacedim>::
3942 cell_relation_t>
3944
3950public:
3961 DeclException2(ExcInvalidLevel,
3962 int,
3963 int,
3964 << "You are requesting information from refinement level "
3965 << arg1
3966 << " of a triangulation, but this triangulation only has "
3967 << arg2 << " refinement levels. The given level " << arg1
3968 << " must be *less* than " << arg2 << '.');
3976 ExcTriangulationNotEmpty,
3977 int,
3978 int,
3979 << "You are trying to perform an operation on a triangulation "
3980 << "that is only allowed if the triangulation is currently empty. "
3981 << "However, it currently stores " << arg1 << " vertices and has "
3982 << "cells on " << arg2 << " levels.");
3988 DeclException0(ExcGridReadError);
3993 DeclException0(ExcFacesHaveNoLevel);
3999 DeclException1(ExcEmptyLevel,
4000 int,
4001 << "You tried to do something on level " << arg1
4002 << ", but this level is empty.");
4003
4011 DeclException1(ExcBoundaryIdNotFound,
4013 << "The given boundary_id " << arg1
4014 << " is not defined in this Triangulation!");
4015
4022 ExcInconsistentCoarseningFlags,
4023 "A cell is flagged for coarsening, but either not all of its siblings "
4024 "are active or flagged for coarsening as well. Please clean up all "
4025 "coarsen flags on your triangulation via "
4026 "Triangulation::prepare_coarsening_and_refinement() beforehand!");
4027
4030protected:
4036
4041 std::vector<ReferenceCell> reference_cells;
4042
4056 static void
4057 write_bool_vector(const unsigned int magic_number1,
4058 const std::vector<bool> &v,
4059 const unsigned int magic_number2,
4060 std::ostream &out);
4061
4066 static void
4067 read_bool_vector(const unsigned int magic_number1,
4068 std::vector<bool> &v,
4069 const unsigned int magic_number2,
4070 std::istream &in);
4071
4076 void
4078
4082 virtual void
4084
4085
4086private:
4091 std::unique_ptr<
4094
4101 std::vector<GridTools::PeriodicFacePair<cell_iterator>>
4103
4108 std::map<std::pair<cell_iterator, unsigned int>,
4109 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>>
4111
4122 TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
4125 using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
4126 using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
4127 using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
4128
4139 begin_raw(const unsigned int level = 0) const;
4140
4146 end_raw(const unsigned int level) const;
4147
4162 begin_raw_line(const unsigned int level = 0) const;
4163
4182 begin_line(const unsigned int level = 0) const;
4183
4202 begin_active_line(const unsigned int level = 0) const;
4203
4209 end_line() const;
4210
4238 begin_raw_quad(const unsigned int level = 0) const;
4239
4258 begin_quad(const unsigned int level = 0) const;
4259
4278 begin_active_quad(const unsigned int level = 0) const;
4279
4285 end_quad() const;
4286
4313 begin_raw_hex(const unsigned int level = 0) const;
4314
4333 begin_hex(const unsigned int level = 0) const;
4334
4353 begin_active_hex(const unsigned int level = 0) const;
4354
4360 end_hex() const;
4361
4378 void
4380
4384 void
4386
4393 void
4395
4399 void
4401
4405 void
4407
4423
4430 void
4432
4437 void
4439
4455 virtual unsigned int
4457 const types::coarse_cell_id coarse_cell_id) const;
4458
4459
4472 virtual types::coarse_cell_id
4474 const unsigned int coarse_cell_index) const;
4475
4480 std::vector<
4481 std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
4483
4489 std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
4491
4492
4496 std::vector<Point<spacedim>> vertices;
4497
4501 std::vector<bool> vertices_used;
4502
4506 std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4508
4513
4514
4520
4531
4546 std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4548
4549
4569 std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4571
4572 // make a couple of classes friends
4573 template <int, int, int>
4574 friend class TriaAccessorBase;
4575 template <int, int, int>
4576 friend class TriaAccessor;
4577 friend class TriaAccessor<0, 1, spacedim>;
4578
4579 friend class CellAccessor<dim, spacedim>;
4580
4581 friend struct ::internal::TriaAccessorImplementation::Implementation;
4582
4583 friend struct ::internal::TriangulationImplementation::Implementation;
4584 friend struct ::internal::TriangulationImplementation::
4585 ImplementationMixedMesh;
4586
4587 friend class ::internal::TriangulationImplementation::TriaObjects;
4588
4589 // explicitly check for sensible template arguments, but not on windows
4590 // because MSVC creates bogus warnings during normal compilation
4591#ifndef DEAL_II_MSVC
4592 static_assert(dim <= spacedim,
4593 "The dimension <dim> of a Triangulation must be less than or "
4594 "equal to the space dimension <spacedim> in which it lives.");
4595#endif
4596};
4597
4598
4599#ifndef DOXYGEN
4600
4601
4602
4603namespace internal
4604{
4605 namespace TriangulationImplementation
4606 {
4607 template <class Archive>
4608 void
4609 NumberCache<1>::serialize(Archive &ar, const unsigned int)
4610 {
4611 ar &n_levels;
4612 ar &n_lines &n_lines_level;
4613 ar &n_active_lines &n_active_lines_level;
4614 }
4615
4616
4617 template <class Archive>
4618 void
4619 NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4620 {
4621 this->NumberCache<1>::serialize(ar, version);
4622
4623 ar &n_quads &n_quads_level;
4624 ar &n_active_quads &n_active_quads_level;
4625 }
4626
4627
4628 template <class Archive>
4629 void
4630 NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4631 {
4632 this->NumberCache<2>::serialize(ar, version);
4633
4634 ar &n_hexes &n_hexes_level;
4635 ar &n_active_hexes &n_active_hexes_level;
4636 }
4637
4638 } // namespace TriangulationImplementation
4639} // namespace internal
4640
4641
4642template <int dim, int spacedim>
4645 const unsigned int index) const
4646{
4647 AssertIndexRange(index, vertices_used.size());
4648 return vertices_used[index];
4649}
4650
4651
4652
4653template <int dim, int spacedim>
4655inline unsigned int Triangulation<dim, spacedim>::n_levels() const
4656{
4657 return number_cache.n_levels;
4658}
4659
4660template <int dim, int spacedim>
4662inline unsigned int Triangulation<dim, spacedim>::n_global_levels() const
4663{
4664 return number_cache.n_levels;
4665}
4666
4667
4668template <int dim, int spacedim>
4670inline unsigned int Triangulation<dim, spacedim>::n_vertices() const
4671{
4672 return vertices.size();
4673}
4674
4675
4676
4677template <int dim, int spacedim>
4679inline const std::vector<Point<spacedim>>
4681{
4682 return vertices;
4683}
4684
4685
4686template <int dim, int spacedim>
4688template <class Archive>
4689void Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4690{
4691 // as discussed in the documentation, do not store the signals as
4692 // well as boundary and manifold description but everything else
4693 ar &smooth_grid;
4694
4695 unsigned int n_levels = levels.size();
4696 ar &n_levels;
4697 for (const auto &level : levels)
4698 ar &level;
4699
4700 // boost dereferences a nullptr when serializing a nullptr
4701 // at least up to 1.65.1. This causes problems with clang-5.
4702 // Therefore, work around it.
4703 bool faces_is_nullptr = (faces.get() == nullptr);
4704 ar &faces_is_nullptr;
4705 if (!faces_is_nullptr)
4706 ar &faces;
4707
4708 ar &vertices;
4709 ar &vertices_used;
4710
4711 ar &anisotropic_refinement;
4712 ar &number_cache;
4713
4714 ar &check_for_distorted_cells;
4715
4716 if (dim == 1)
4717 {
4718 ar &vertex_to_boundary_id_map_1d;
4719 ar &vertex_to_manifold_id_map_1d;
4720 }
4721}
4722
4723
4724
4725template <int dim, int spacedim>
4727template <class Archive>
4728void Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4729{
4730 // clear previous content. this also calls the respective signal
4731 clear();
4732
4733 // as discussed in the documentation, do not store the signals as
4734 // well as boundary and manifold description but everything else
4735 ar &smooth_grid;
4736
4737 unsigned int size;
4738 ar &size;
4739 levels.resize(size);
4740 for (auto &level_ : levels)
4741 {
4742 std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4743 ar &level;
4744 level_ = std::move(level);
4745 }
4746
4747 // Workaround for nullptr, see in save().
4748 bool faces_is_nullptr = true;
4749 ar &faces_is_nullptr;
4750 if (!faces_is_nullptr)
4751 ar &faces;
4752
4753 ar &vertices;
4754 ar &vertices_used;
4755
4756 ar &anisotropic_refinement;
4757 ar &number_cache;
4758
4759 // the levels do not serialize the active_cell_indices because
4760 // they are easy enough to rebuild upon re-loading data. do
4761 // this here. don't forget to first resize the fields appropriately
4762 {
4763 for (const auto &level : levels)
4764 {
4765 level->active_cell_indices.resize(level->refine_flags.size());
4766 level->global_active_cell_indices.resize(level->refine_flags.size());
4767 level->global_level_cell_indices.resize(level->refine_flags.size());
4768 }
4769 reset_cell_vertex_indices_cache();
4770 reset_active_cell_indices();
4771 reset_global_cell_indices();
4772 }
4773
4774 reset_policy();
4775
4776 bool my_check_for_distorted_cells;
4777 ar &my_check_for_distorted_cells;
4778
4779 Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4780 ExcMessage("The triangulation loaded into here must have the "
4781 "same setting with regard to reporting distorted "
4782 "cell as the one previously stored."));
4783
4784 if (dim == 1)
4785 {
4786 ar &vertex_to_boundary_id_map_1d;
4787 ar &vertex_to_manifold_id_map_1d;
4788 }
4789
4790 // trigger the create signal to indicate
4791 // that new content has been imported into
4792 // the triangulation
4793 signals.create();
4794}
4795
4796
4797
4798template <int dim, int spacedim>
4800inline unsigned int Triangulation<dim, spacedim>::
4802 const types::coarse_cell_id coarse_cell_id) const
4803{
4804 return coarse_cell_id;
4805}
4806
4807
4808
4809template <int dim, int spacedim>
4813 const unsigned int coarse_cell_index) const
4814{
4815 return coarse_cell_index;
4816}
4817
4818
4819
4820/* -------------- declaration of explicit specializations ------------- */
4821
4822template <>
4823unsigned int
4825template <>
4826unsigned int
4827Triangulation<1, 1>::n_quads(const unsigned int level) const;
4828template <>
4829unsigned int
4830Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4831template <>
4832unsigned int
4833Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4834template <>
4835unsigned int
4836Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4837template <>
4838unsigned int
4840template <>
4841unsigned int
4842Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4843template <>
4844unsigned int
4846template <>
4847unsigned int
4848Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4849template <>
4850unsigned int
4851Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4852template <>
4853unsigned int
4855template <>
4856unsigned int
4858template <>
4859unsigned int
4860Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4861template <>
4862unsigned int
4863Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4864
4865template <>
4866unsigned int
4868
4869
4870// -------------------------------------------------------------------
4871// -- Explicit specializations for codimension one grids
4872
4873
4874template <>
4875unsigned int
4877template <>
4878unsigned int
4879Triangulation<1, 2>::n_quads(const unsigned int level) const;
4880template <>
4881unsigned int
4882Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4883template <>
4884unsigned int
4885Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4886template <>
4887unsigned int
4888Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4889template <>
4890unsigned int
4891Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4892template <>
4893unsigned int
4895template <>
4896unsigned int
4898
4899// -------------------------------------------------------------------
4900// -- Explicit specializations for codimension two grids
4901
4902template <>
4903unsigned int
4905template <>
4906unsigned int
4907Triangulation<1, 3>::n_quads(const unsigned int level) const;
4908template <>
4909unsigned int
4910Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4911template <>
4912unsigned int
4913Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4914template <>
4915unsigned int
4916Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4917template <>
4918unsigned int
4919Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4920template <>
4921unsigned int
4923template <>
4924unsigned int
4926
4927template <>
4928bool
4930template <>
4931bool
4933template <>
4934bool
4936
4937
4938extern template class Triangulation<1, 1>;
4939extern template class Triangulation<1, 2>;
4940extern template class Triangulation<1, 3>;
4941extern template class Triangulation<2, 2>;
4942extern template class Triangulation<2, 3>;
4943extern template class Triangulation<3, 3>;
4944
4945#endif // DOXYGEN
4946
4948
4949// Include tria_accessor.h here, so that it is possible for an end
4950// user to use the iterators of Triangulation<dim> directly without
4951// the need to include tria_accessor.h separately. (Otherwise the
4952// iterators are an 'opaque' or 'incomplete' type.)
4954
4955#endif
CellStatus
Definition cell_status.h:31
@ cell_will_be_refined
@ children_will_be_coarsened
Definition point.h:111
unsigned int n_hexs(const unsigned int level) const
virtual void add_periodicity(const std::vector< GridTools::PeriodicFacePair< cell_iterator > > &)
void save_user_flags_line(std::vector< bool > &v) const
virtual types::global_cell_index n_global_active_cells() const
quad_iterator begin_quad(const unsigned int level=0) const
typename IteratorSelector::raw_line_iterator raw_line_iterator
Definition tria.h:4125
active_vertex_iterator begin_active_vertex() const
void load_user_indices_quad(const std::vector< unsigned int > &v)
unsigned int n_quads() const
void load_user_indices(const std::vector< unsigned int > &v)
std::vector< bool > vertices_used
Definition tria.h:4501
bool anisotropic_refinement
Definition tria.h:4512
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, unsigned char > > periodic_face_map
Definition tria.h:4110
active_quad_iterator begin_active_quad(const unsigned int level=0) const
bool get_anisotropic_refinement_flag() const
virtual types::coarse_cell_id n_global_coarse_cells() const
std::unique_ptr< std::map< unsigned int, types::manifold_id > > vertex_to_manifold_id_map_1d
Definition tria.h:4570
void save_user_pointers_quad(std::vector< void * > &v) const
void save_user_flags_hex(std::ostream &out) const
void clear_user_flags_quad()
unsigned int n_faces() const
active_hex_iterator begin_active_hex(const unsigned int level=0) const
static void read_bool_vector(const unsigned int magic_number1, std::vector< bool > &v, const unsigned int magic_number2, std::istream &in)
unsigned int n_active_lines(const unsigned int level) const
bool all_reference_cells_are_hyper_cube() const
void load_user_flags_line(std::istream &in)
void clear_user_data()
raw_hex_iterator begin_raw_hex(const unsigned int level=0) const
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, unsigned char > > & get_periodic_face_map() const
void save_user_flags_line(std::ostream &out) const
active_cell_iterator last_active() const
void save(Archive &ar, const unsigned int version) const
std::vector< typename internal::CellAttachedDataSerializer< dim, spacedim >::cell_relation_t > local_cell_relations
Definition tria.h:3943
const Triangulation< dim, spacedim > & get_triangulation() const
void reset_global_cell_indices()
virtual void save(const std::string &file_basename) const
face_iterator end_face() const
void reset_active_cell_indices()
cell_iterator create_cell_iterator(const CellId &cell_id) const
cell_iterator begin(const unsigned int level=0) const
void fix_coarsen_flags()
virtual void load(const std::string &file_basename)
unsigned int n_active_cells(const unsigned int level) const
unsigned int n_cells(const unsigned int level) const
void save_user_pointers_line(std::vector< void * > &v) const
void save_user_flags_quad(std::vector< bool > &v) const
void load_refine_flags(std::istream &in)
void save_user_indices_line(std::vector< unsigned int > &v) const
raw_cell_iterator begin_raw(const unsigned int level=0) const
unsigned int n_lines() const
void load_coarsen_flags(const std::vector< bool > &v)
unsigned int n_raw_lines() const
virtual std::size_t memory_consumption() const
std::vector< Point< spacedim > > vertices
Definition tria.h:4496
raw_quad_iterator begin_raw_quad(const unsigned int level=0) const
virtual types::subdomain_id locally_owned_subdomain() const
unsigned int n_raw_faces() const
unsigned int n_active_faces() const
const bool check_for_distorted_cells
Definition tria.h:4519
raw_cell_iterator end_raw(const unsigned int level) const
line_iterator end_line() const
std::unique_ptr< std::map< unsigned int, types::boundary_id > > vertex_to_boundary_id_map_1d
Definition tria.h:4547
void load_user_flags_quad(std::istream &in)
unsigned int n_active_cells() const
virtual void update_reference_cells()
std::vector< ReferenceCell > reference_cells
Definition tria.h:4041
void update_periodic_face_map()
unsigned int n_raw_lines(const unsigned int level) const
void clear_despite_subscriptions()
void save_user_flags(std::ostream &out) const
void load_user_flags_hex(std::istream &in)
const std::vector< Point< spacedim > > & get_vertices() const
void load_user_pointers_quad(const std::vector< void * > &v)
std::unique_ptr<::internal::TriangulationImplementation::TriaFaces > faces
Definition tria.h:4490
::CellStatus CellStatus
Definition tria.h:2241
unsigned int n_used_vertices() const
void reset_cell_vertex_indices_cache()
unsigned int n_active_lines() const
unsigned int n_levels() const
void load_user_flags_hex(const std::vector< bool > &v)
void load_user_indices_line(const std::vector< unsigned int > &v)
void clear_user_flags_hex()
void save_user_pointers_hex(std::vector< void * > &v) const
const std::vector< ReferenceCell > & get_reference_cells() const
typename IteratorSelector::raw_quad_iterator raw_quad_iterator
Definition tria.h:4126
void load_user_pointers(const std::vector< void * > &v)
unsigned int register_data_attach(const std::function< std::vector< char >(const cell_iterator &, const ::CellStatus)> &pack_callback, const bool returns_variable_size_data)
::internal::TriangulationImplementation::NumberCache< dim > number_cache
Definition tria.h:4530
void load_user_flags_quad(const std::vector< bool > &v)
void save_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &file_basename) const
void load_user_flags_line(const std::vector< bool > &v)
void save_user_indices_hex(std::vector< unsigned int > &v) const
DistortedCellList execute_refinement()
void update_cell_relations()
active_line_iterator begin_active_line(const unsigned int level=0) const
void save_user_indices_quad(std::vector< unsigned int > &v) const
void load_user_pointers_hex(const std::vector< void * > &v)
void pack_data_serial()
cell_iterator end() const
virtual bool has_hanging_nodes() const
std::vector< GridTools::PeriodicFacePair< cell_iterator > > periodic_face_pairs_level_0
Definition tria.h:4102
unsigned int n_raw_cells(const unsigned int level) const
cell_iterator end(const unsigned int level) const
bool contains_cell(const CellId &cell_id) const
void load_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const unsigned int local_num_cells, const std::string &file_basename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable)
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const
unsigned int n_raw_quads(const unsigned int level) const
internal::CellAttachedData< dim, spacedim > cell_attached_data
Definition tria.h:3877
void load_coarsen_flags(std::istream &out)
quad_iterator end_quad() const
line_iterator begin_line(const unsigned int level=0) const
unsigned int max_adjacent_cells() const
vertex_iterator begin_vertex() const
void clear_user_flags()
unsigned int n_hexs() const
vertex_iterator end_vertex() const
bool vertex_used(const unsigned int index) const
void load_user_pointers_line(const std::vector< void * > &v)
void save_user_flags(std::vector< bool > &v) const
hex_iterator end_hex() const
hex_iterator begin_hex(const unsigned int level=0) const
virtual unsigned int n_global_levels() const
active_cell_iterator end_active(const unsigned int level) const
bool is_mixed_mesh() const
cell_iterator last() const
unsigned int n_active_quads() const
void save_coarsen_flags(std::vector< bool > &v) const
void load_user_indices_hex(const std::vector< unsigned int > &v)
unsigned int n_raw_quads() const
void save_user_pointers(std::vector< void * > &v) const
face_iterator begin_face() const
unsigned int n_cells() const
virtual bool prepare_coarsening_and_refinement()
unsigned int n_active_quads(const unsigned int level) const
void unpack_data_serial()
const std::vector< bool > & get_used_vertices() const
void load_user_flags(const std::vector< bool > &v)
typename IteratorSelector::raw_hex_iterator raw_hex_iterator
Definition tria.h:4127
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
Definition tria.h:4507
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition tria.h:4035
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition tria.h:4093
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition tria.h:2526
internal::CellAttachedDataSerializer< dim, spacedim > data_serializer
Definition tria.h:3945
unsigned int n_active_hexs(const unsigned int level) const
unsigned int n_vertices() const
void load(Archive &ar, const unsigned int version)
void save_user_indices(std::vector< unsigned int > &v) const
void notify_ready_to_unpack(const unsigned int handle, const std::function< void(const cell_iterator &, const ::CellStatus, const boost::iterator_range< std::vector< char >::const_iterator > &)> &unpack_callback)
void save_user_flags_hex(std::vector< bool > &v) const
bool all_reference_cells_are_simplex() const
std::vector< std::unique_ptr<::internal::TriangulationImplementation::TriaLevel > > levels
Definition tria.h:4482
unsigned int n_raw_hexs(const unsigned int level) const
unsigned int n_lines(const unsigned int level) const
unsigned int n_active_hexs() const
virtual unsigned int coarse_cell_id_to_coarse_cell_index(const types::coarse_cell_id coarse_cell_id) const
unsigned int n_quads(const unsigned int level) const
void load_user_flags(std::istream &in)
void reset_policy()
void save_coarsen_flags(std::ostream &out) const
active_face_iterator begin_active_face() const
void clear_user_flags_line()
raw_line_iterator begin_raw_line(const unsigned int level=0) const
static void write_bool_vector(const unsigned int magic_number1, const std::vector< bool > &v, const unsigned int magic_number2, std::ostream &out)
active_cell_iterator begin_active(const unsigned int level=0) const
void execute_coarsening()
void save_refine_flags(std::vector< bool > &v) const
std::vector< char > dest_data_variable
Definition tria.h:535
std::vector< char > src_data_fixed
Definition tria.h:525
std::vector< char > dest_data_fixed
Definition tria.h:526
std::vector< int > src_sizes_variable
Definition tria.h:532
std::vector< unsigned int > sizes_fixed_cumulative
Definition tria.h:519
typename std::pair< cell_iterator, CellStatus > cell_relation_t
Definition tria.h:393
std::vector< int > dest_sizes_variable
Definition tria.h:533
std::vector< char > src_data_variable
Definition tria.h:534
#define DEAL_II_DEPRECATED
Definition config.h:205
#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
Point< 2 > first
Definition grid_out.cc:4629
unsigned int level
Definition grid_out.cc:4632
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
IteratorRange< active_face_iterator > active_face_iterators() const
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
IteratorRange< cell_iterator > cell_iterators() const
#define DeclException0(Exception0)
Definition exceptions.h:466
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:534
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:489
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:511
static ::ExceptionBase & ExcMessage(std::string arg1)
typename IteratorSelector::hex_iterator hex_iterator
Definition tria.h:1691
typename IteratorSelector::active_quad_iterator active_quad_iterator
Definition tria.h:1682
typename IteratorSelector::active_hex_iterator active_hex_iterator
Definition tria.h:1702
typename IteratorSelector::quad_iterator quad_iterator
Definition tria.h:1667
typename IteratorSelector::line_iterator line_iterator
Definition tria.h:1643
typename IteratorSelector::active_line_iterator active_line_iterator
Definition tria.h:1658
std::size_t size
Definition mpi.cc:734
STL namespace.
Definition types.h:32
global_cell_index coarse_cell_id
Definition types.h:130
T operator()(InputIterator first, InputIterator last) const
Definition tria.h:2284
virtual ~DistortedCellList() noexcept override
boost::signals2::signal< void()> mesh_movement
Definition tria.h:2341
boost::signals2::signal< void()> post_distributed_load
Definition tria.h:2515
boost::signals2::signal< void()> pre_distributed_save
Definition tria.h:2493
boost::signals2::signal< unsigned int(const cell_iterator &, const ::CellStatus), CellWeightSum< unsigned int > > weight
Definition tria.h:2439
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition tria.h:2372
boost::signals2::signal< void()> pre_distributed_refinement
Definition tria.h:2452
boost::signals2::signal< void()> post_distributed_refinement
Definition tria.h:2471
boost::signals2::signal< void()> pre_refinement
Definition tria.h:2317
boost::signals2::signal< void()> any_change
Definition tria.h:2398
boost::signals2::signal< void()> create
Definition tria.h:2308
boost::signals2::signal< void()> clear
Definition tria.h:2387
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition tria.h:2362
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition tria.h:2352
boost::signals2::signal< void()> post_refinement
Definition tria.h:2324
boost::signals2::signal< void()> pre_partition
Definition tria.h:2332
boost::signals2::signal< void()> pre_distributed_repartition
Definition tria.h:2478
boost::signals2::signal< void()> pre_distributed_load
Definition tria.h:2508
boost::signals2::signal< void()> post_p4est_refinement
Definition tria.h:2461
boost::signals2::signal< void()> post_distributed_repartition
Definition tria.h:2485
boost::signals2::signal< void()> post_distributed_save
Definition tria.h:2500
std::vector< pack_callback_t > pack_callbacks_variable
Definition tria.h:362
unsigned int n_attached_deserialize
Definition tria.h:352
std::vector< pack_callback_t > pack_callbacks_fixed
Definition tria.h:361
unsigned int n_attached_data_sets
Definition tria.h:346
std::function< std::vector< char >(cell_iterator, CellStatus)> pack_callback_t
Definition tria.h:355
std::shared_ptr< const Utilities::MPI::Partitioner > active_cell_index_partitioner
Definition tria.h:184
void serialize(Archive &ar, const unsigned int version)
std::vector< unsigned int > n_active_lines_level
Definition tria.h:178
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > level_cell_index_partitioners
Definition tria.h:190
std::vector< unsigned int > n_active_quads_level
Definition tria.h:248
void serialize(Archive &ar, const unsigned int version)
void serialize(Archive &ar, const unsigned int version)
std::vector< unsigned int > n_active_hexes_level
Definition tria.h:307