Reference documentation for deal.II version GIT relicensing-426-g7976cfd195 2024-04-18 21:10:01+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 - 2023 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
24#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
387 using cell_relation_t = typename std::pair<cell_iterator, CellStatus>;
388
390
399 void
400 pack_data(
401 const std::vector<cell_relation_t> &cell_relations,
402 const std::vector<
404 &pack_callbacks_fixed,
405 const std::vector<
407 &pack_callbacks_variable,
408 const MPI_Comm &mpi_communicator);
409
417 void
418 unpack_cell_status(std::vector<cell_relation_t> &cell_relations) const;
419
432 void
433 unpack_data(
434 const std::vector<cell_relation_t> &cell_relations,
435 const unsigned int handle,
436 const std::function<
437 void(const cell_iterator &,
438 const CellStatus &,
439 const boost::iterator_range<std::vector<char>::const_iterator> &)>
440 &unpack_callback) const;
441
456 void
457 save(const unsigned int global_first_cell,
458 const unsigned int global_num_cells,
459 const std::string &filename,
460 const MPI_Comm &mpi_communicator) const;
461
480 void
481 load(const unsigned int global_first_cell,
482 const unsigned int global_num_cells,
483 const unsigned int local_num_cells,
484 const std::string &filename,
485 const unsigned int n_attached_deserialize_fixed,
486 const unsigned int n_attached_deserialize_variable,
487 const MPI_Comm &mpi_communicator);
488
495 void
496 clear();
497
502
513 std::vector<unsigned int> sizes_fixed_cumulative;
514
519 std::vector<char> src_data_fixed;
520 std::vector<char> dest_data_fixed;
521
526 std::vector<int> src_sizes_variable;
527 std::vector<int> dest_sizes_variable;
528 std::vector<char> src_data_variable;
529 std::vector<char> dest_data_variable;
530 };
531} // namespace internal
532
533
534/*------------------------------------------------------------------------*/
535
536
1314template <int dim, int spacedim = dim>
1317{
1318private:
1325
1326public:
1332 {
1337 none = 0x0,
1381 limit_level_difference_at_vertices = 0x1,
1402 eliminate_unrefined_islands = 0x2,
1418 patch_level_1 = 0x4,
1442 coarsest_level_1 = 0x8,
1467 allow_anisotropic_smoothing = 0x10,
1500 eliminate_refined_inner_islands = 0x100,
1505 eliminate_refined_boundary_islands = 0x200,
1511 do_not_produce_unrefined_islands = 0x400,
1512
1517 smoothing_on_refinement =
1518 (limit_level_difference_at_vertices | eliminate_unrefined_islands),
1523 smoothing_on_coarsening =
1524 (eliminate_refined_inner_islands | eliminate_refined_boundary_islands |
1525 do_not_produce_unrefined_islands),
1526
1532 maximum_smoothing = 0xffff ^ allow_anisotropic_smoothing
1534
1551
1557
1575
1589 using face_iterator = TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1590
1603 TriaActiveIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1604
1614
1629
1637 using line_iterator = typename IteratorSelector::line_iterator;
1638
1652 using active_line_iterator = typename IteratorSelector::active_line_iterator;
1653
1661 using quad_iterator = typename IteratorSelector::quad_iterator;
1662
1676 using active_quad_iterator = typename IteratorSelector::active_quad_iterator;
1677
1685 using hex_iterator = typename IteratorSelector::hex_iterator;
1686
1696 using active_hex_iterator = typename IteratorSelector::active_hex_iterator;
1697
1717 {
1724 virtual ~DistortedCellList() noexcept override;
1725
1730 std::list<typename Triangulation<dim, spacedim>::cell_iterator>
1731 distorted_cells;
1732 };
1733
1737 static constexpr unsigned int dimension = dim;
1738
1742 static constexpr unsigned int space_dimension = spacedim;
1743
1758 Triangulation(const MeshSmoothing smooth_grid = none,
1759 const bool check_for_distorted_cells = false);
1760
1776 Triangulation(const Triangulation<dim, spacedim> &) = delete;
1777
1784 Triangulation(Triangulation<dim, spacedim> &&tria) noexcept;
1785
1790 operator=(Triangulation<dim, spacedim> &&tria) noexcept;
1791
1795 virtual ~Triangulation() override;
1796
1803 virtual void
1804 clear();
1805
1810 virtual MPI_Comm
1811 get_communicator() const;
1812
1818 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1819 global_active_cell_index_partitioner() const;
1820
1826 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1827 global_level_cell_index_partitioner(const unsigned int level) const;
1828
1833 virtual void
1834 set_mesh_smoothing(const MeshSmoothing mesh_smoothing);
1835
1839 virtual const MeshSmoothing &
1840 get_mesh_smoothing() const;
1841
1864 void
1865 set_manifold(const types::manifold_id number,
1866 const Manifold<dim, spacedim> &manifold_object);
1867
1884 void
1885 reset_manifold(const types::manifold_id manifold_number);
1886
1903 void
1904 reset_all_manifolds();
1905
1914 void
1915 set_all_manifold_ids(const types::manifold_id number);
1916
1925 void
1926 set_all_manifold_ids_on_boundary(const types::manifold_id number);
1927
1937 void
1938 set_all_manifold_ids_on_boundary(const types::boundary_id b_id,
1939 const types::manifold_id number);
1940
1952 const Manifold<dim, spacedim> &
1953 get_manifold(const types::manifold_id number) const;
1954
1967 virtual std::vector<types::boundary_id>
1968 get_boundary_ids() const;
1969
1982 virtual std::vector<types::manifold_id>
1983 get_manifold_ids() const;
1984
2006 virtual void
2007 copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
2008
2057 virtual void
2058 create_triangulation(const std::vector<Point<spacedim>> &vertices,
2059 const std::vector<CellData<dim>> &cells,
2060 const SubCellData &subcelldata);
2061
2074 virtual void
2075 create_triangulation(
2076 const TriangulationDescription::Description<dim, spacedim>
2077 &construction_data);
2078
2087 void
2088 flip_all_direction_flags();
2089
2101 void
2102 set_all_refine_flags();
2103
2126 void
2127 refine_global(const unsigned int times = 1);
2128
2141 void
2142 coarsen_global(const unsigned int times = 1);
2143
2175 virtual void
2176 execute_coarsening_and_refinement();
2177
2207 virtual bool
2208 prepare_coarsening_and_refinement();
2209
2224 using CellStatus DEAL_II_DEPRECATED_EARLY = ::CellStatus;
2225
2230 static constexpr auto CELL_PERSIST DEAL_II_DEPRECATED_EARLY =
2232
2237 static constexpr auto CELL_REFINE DEAL_II_DEPRECATED_EARLY =
2239
2244 static constexpr auto CELL_COARSEN DEAL_II_DEPRECATED_EARLY =
2246
2251 static constexpr auto CELL_INVALID DEAL_II_DEPRECATED_EARLY =
2253
2254
2260 template <typename T>
2262 {
2263 using result_type = T;
2264
2265 template <typename InputIterator>
2266 T
2267 operator()(InputIterator first, InputIterator last) const
2268 {
2269 return std::accumulate(first, last, T());
2270 }
2271 };
2272
2282 struct Signals
2283 {
2291 boost::signals2::signal<void()> create;
2292
2300 boost::signals2::signal<void()> pre_refinement;
2301
2307 boost::signals2::signal<void()> post_refinement;
2308
2315 boost::signals2::signal<void()> pre_partition;
2316
2324 boost::signals2::signal<void()> mesh_movement;
2325
2333 boost::signals2::signal<void(
2334 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2336
2343 boost::signals2::signal<void(
2344 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2346
2353 boost::signals2::signal<void(
2354 const Triangulation<dim, spacedim> &destination_tria)>
2356
2370 boost::signals2::signal<void()> clear;
2371
2381 boost::signals2::signal<void()> any_change;
2382
2419 boost::signals2::signal<unsigned int(const cell_iterator &,
2420 const ::CellStatus),
2423
2435 boost::signals2::signal<void()> pre_distributed_refinement;
2436
2444 boost::signals2::signal<void()> post_p4est_refinement;
2445
2454 boost::signals2::signal<void()> post_distributed_refinement;
2455
2461 boost::signals2::signal<void()> pre_distributed_repartition;
2462
2468 boost::signals2::signal<void()> post_distributed_repartition;
2469
2476 boost::signals2::signal<void()> pre_distributed_save;
2477
2483 boost::signals2::signal<void()> post_distributed_save;
2484
2491 boost::signals2::signal<void()> pre_distributed_load;
2492
2498 boost::signals2::signal<void()> post_distributed_load;
2499 };
2500
2510
2522 void
2523 save_refine_flags(std::ostream &out) const;
2524
2528 void
2529 save_refine_flags(std::vector<bool> &v) const;
2530
2534 void
2535 load_refine_flags(std::istream &in);
2536
2540 void
2541 load_refine_flags(const std::vector<bool> &v);
2542
2546 void
2547 save_coarsen_flags(std::ostream &out) const;
2548
2552 void
2553 save_coarsen_flags(std::vector<bool> &v) const;
2554
2558 void
2559 load_coarsen_flags(std::istream &out);
2560
2564 void
2565 load_coarsen_flags(const std::vector<bool> &v);
2566
2571 bool
2573
2585 void
2587
2593 void
2594 save_user_flags(std::ostream &out) const;
2595
2601 void
2602 save_user_flags(std::vector<bool> &v) const;
2603
2608 void
2609 load_user_flags(std::istream &in);
2610
2615 void
2616 load_user_flags(const std::vector<bool> &v);
2617
2622 void
2624
2629 void
2630 save_user_flags_line(std::ostream &out) const;
2631
2637 void
2638 save_user_flags_line(std::vector<bool> &v) const;
2639
2644 void
2645 load_user_flags_line(std::istream &in);
2646
2651 void
2652 load_user_flags_line(const std::vector<bool> &v);
2653
2658 void
2660
2665 void
2666 save_user_flags_quad(std::ostream &out) const;
2667
2673 void
2674 save_user_flags_quad(std::vector<bool> &v) const;
2675
2680 void
2681 load_user_flags_quad(std::istream &in);
2682
2687 void
2688 load_user_flags_quad(const std::vector<bool> &v);
2689
2690
2695 void
2697
2702 void
2703 save_user_flags_hex(std::ostream &out) const;
2704
2710 void
2711 save_user_flags_hex(std::vector<bool> &v) const;
2712
2717 void
2718 load_user_flags_hex(std::istream &in);
2719
2724 void
2725 load_user_flags_hex(const std::vector<bool> &v);
2726
2732 void
2734
2740 void
2741 save_user_indices(std::vector<unsigned int> &v) const;
2742
2747 void
2748 load_user_indices(const std::vector<unsigned int> &v);
2749
2755 void
2756 save_user_pointers(std::vector<void *> &v) const;
2757
2762 void
2763 load_user_pointers(const std::vector<void *> &v);
2764
2770 void
2771 save_user_indices_line(std::vector<unsigned int> &v) const;
2772
2777 void
2778 load_user_indices_line(const std::vector<unsigned int> &v);
2779
2785 void
2786 save_user_indices_quad(std::vector<unsigned int> &v) const;
2787
2792 void
2793 load_user_indices_quad(const std::vector<unsigned int> &v);
2794
2800 void
2801 save_user_indices_hex(std::vector<unsigned int> &v) const;
2802
2807 void
2808 load_user_indices_hex(const std::vector<unsigned int> &v);
2814 void
2815 save_user_pointers_line(std::vector<void *> &v) const;
2816
2821 void
2822 load_user_pointers_line(const std::vector<void *> &v);
2823
2829 void
2830 save_user_pointers_quad(std::vector<void *> &v) const;
2831
2836 void
2837 load_user_pointers_quad(const std::vector<void *> &v);
2838
2844 void
2845 save_user_pointers_hex(std::vector<void *> &v) const;
2846
2851 void
2852 load_user_pointers_hex(const std::vector<void *> &v);
2853
2879 begin(const unsigned int level = 0) const;
2880
2912 begin_active(const unsigned int level = 0) const;
2913
2919 end() const;
2920
2940 end(const unsigned int level) const;
2941
2962 end_active(const unsigned int level) const;
2963
2964
2969 last() const;
2970
2976
2985 create_cell_iterator(const CellId &cell_id) const;
2986
2997 bool
2998 contains_cell(const CellId &cell_id) const;
3018
3056
3073 cell_iterators_on_level(const unsigned int level) const;
3074
3091 active_cell_iterators_on_level(const unsigned int level) const;
3092
3095 /*-------------------------------------------------------------------------*/
3096
3106 begin_face() const;
3107
3113
3119 end_face() const;
3120
3141
3144 /*-------------------------------------------------------------------------*/
3145
3157
3165
3172 end_vertex() const;
3173
3193 unsigned int
3194 n_lines() const;
3195
3199 unsigned int
3200 n_lines(const unsigned int level) const;
3201
3205 unsigned int
3207
3211 unsigned int
3212 n_active_lines(const unsigned int level) const;
3213
3217 unsigned int
3218 n_quads() const;
3219
3223 unsigned int
3224 n_quads(const unsigned int level) const;
3225
3229 unsigned int
3231
3235 unsigned int
3236 n_active_quads(const unsigned int level) const;
3237
3241 unsigned int
3242 n_hexs() const;
3243
3248 unsigned int
3249 n_hexs(const unsigned int level) const;
3250
3254 unsigned int
3256
3261 unsigned int
3262 n_active_hexs(const unsigned int level) const;
3263
3268 unsigned int
3269 n_cells() const;
3270
3275 unsigned int
3276 n_cells(const unsigned int level) const;
3277
3282 unsigned int
3284
3294
3295
3300 unsigned int
3301 n_active_cells(const unsigned int level) const;
3302
3307 virtual types::coarse_cell_id
3309
3315 unsigned int
3316 n_faces() const;
3317
3323 unsigned int
3325
3343 unsigned int
3344 n_levels() const;
3345
3352 virtual unsigned int
3354
3364 virtual bool
3366
3374 unsigned int
3375 n_vertices() const;
3376
3385 const std::vector<Point<spacedim>> &
3387
3392 unsigned int
3394
3398 bool
3399 vertex_used(const unsigned int index) const;
3400
3405 const std::vector<bool> &
3407
3419 unsigned int
3421
3428 virtual types::subdomain_id
3430
3442
3449
3450
3467 unsigned int
3469
3479 unsigned int
3480 n_raw_lines(const unsigned int level) const;
3481
3491 unsigned int
3493
3503 unsigned int
3504 n_raw_quads(const unsigned int level) const;
3505
3515 unsigned int
3516 n_raw_hexs(const unsigned int level) const;
3517
3527 unsigned int
3528 n_raw_cells(const unsigned int level) const;
3529
3541 unsigned int
3543
3554 virtual std::size_t
3556
3566 template <class Archive>
3567 void
3568 save(Archive &ar, const unsigned int version) const;
3569
3587 template <class Archive>
3588 void
3589 load(Archive &ar, const unsigned int version);
3590
3591
3598 virtual void
3599 save(const std::string &filename) const;
3600
3604 virtual void
3605 load(const std::string &filename);
3606
3607
3623 virtual void
3625 const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3626
3630 const std::map<
3631 std::pair<cell_iterator, unsigned int>,
3632 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>> &
3634
3639 const std::vector<ReferenceCell> &
3641
3646 bool
3648
3653 bool
3655
3661 bool
3663
3664#ifdef DOXYGEN
3670 template <class Archive>
3671 void
3672 serialize(Archive &archive, const unsigned int version);
3673#else
3674 // This macro defines the serialize() method that is compatible with
3675 // the templated save() and load() method that have been implemented.
3676 BOOST_SERIALIZATION_SPLIT_MEMBER()
3677#endif
3678
3683public:
3791 unsigned int
3793 const std::function<std::vector<char>(const cell_iterator &,
3794 const ::CellStatus)>
3795 &pack_callback,
3796 const bool returns_variable_size_data);
3797
3847 void
3849 const unsigned int handle,
3850 const std::function<
3851 void(const cell_iterator &,
3852 const ::CellStatus,
3853 const boost::iterator_range<std::vector<char>::const_iterator> &)>
3854 &unpack_callback);
3855
3857
3858protected:
3865 void
3866 save_attached_data(const unsigned int global_first_cell,
3867 const unsigned int global_num_cells,
3868 const std::string &filename) const;
3869
3877 void
3878 load_attached_data(const unsigned int global_first_cell,
3879 const unsigned int global_num_cells,
3880 const unsigned int local_num_cells,
3881 const std::string &filename,
3882 const unsigned int n_attached_deserialize_fixed,
3883 const unsigned int n_attached_deserialize_variable);
3884
3895 virtual void
3898
3904 std::vector<typename internal::CellAttachedDataSerializer<dim, spacedim>::
3905 cell_relation_t>
3907
3913public:
3924 DeclException2(ExcInvalidLevel,
3925 int,
3926 int,
3927 << "You are requesting information from refinement level "
3928 << arg1
3929 << " of a triangulation, but this triangulation only has "
3930 << arg2 << " refinement levels. The given level " << arg1
3931 << " must be *less* than " << arg2 << '.');
3939 ExcTriangulationNotEmpty,
3940 int,
3941 int,
3942 << "You are trying to perform an operation on a triangulation "
3943 << "that is only allowed if the triangulation is currently empty. "
3944 << "However, it currently stores " << arg1 << " vertices and has "
3945 << "cells on " << arg2 << " levels.");
3951 DeclException0(ExcGridReadError);
3956 DeclException0(ExcFacesHaveNoLevel);
3962 DeclException1(ExcEmptyLevel,
3963 int,
3964 << "You tried to do something on level " << arg1
3965 << ", but this level is empty.");
3966
3974 DeclException1(ExcBoundaryIdNotFound,
3976 << "The given boundary_id " << arg1
3977 << " is not defined in this Triangulation!");
3978
3985 ExcInconsistentCoarseningFlags,
3986 "A cell is flagged for coarsening, but either not all of its siblings "
3987 "are active or flagged for coarsening as well. Please clean up all "
3988 "coarsen flags on your triangulation via "
3989 "Triangulation::prepare_coarsening_and_refinement() beforehand!");
3990
3993protected:
3999
4004 std::vector<ReferenceCell> reference_cells;
4005
4019 static void
4020 write_bool_vector(const unsigned int magic_number1,
4021 const std::vector<bool> &v,
4022 const unsigned int magic_number2,
4023 std::ostream &out);
4024
4029 static void
4030 read_bool_vector(const unsigned int magic_number1,
4031 std::vector<bool> &v,
4032 const unsigned int magic_number2,
4033 std::istream &in);
4034
4039 void
4041
4045 virtual void
4047
4048
4049private:
4054 std::unique_ptr<
4057
4064 std::vector<GridTools::PeriodicFacePair<cell_iterator>>
4066
4071 std::map<std::pair<cell_iterator, unsigned int>,
4072 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>>
4074
4085 TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
4088 using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
4089 using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
4090 using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
4091
4102 begin_raw(const unsigned int level = 0) const;
4103
4109 end_raw(const unsigned int level) const;
4110
4125 begin_raw_line(const unsigned int level = 0) const;
4126
4145 begin_line(const unsigned int level = 0) const;
4146
4165 begin_active_line(const unsigned int level = 0) const;
4166
4172 end_line() const;
4173
4201 begin_raw_quad(const unsigned int level = 0) const;
4202
4221 begin_quad(const unsigned int level = 0) const;
4222
4241 begin_active_quad(const unsigned int level = 0) const;
4242
4248 end_quad() const;
4249
4276 begin_raw_hex(const unsigned int level = 0) const;
4277
4296 begin_hex(const unsigned int level = 0) const;
4297
4316 begin_active_hex(const unsigned int level = 0) const;
4317
4323 end_hex() const;
4324
4341 void
4343
4347 void
4349
4356 void
4358
4362 void
4364
4368 void
4370
4386
4393 void
4395
4400 void
4402
4418 virtual unsigned int
4420 const types::coarse_cell_id coarse_cell_id) const;
4421
4422
4435 virtual types::coarse_cell_id
4437 const unsigned int coarse_cell_index) const;
4438
4443 std::vector<
4444 std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
4446
4452 std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
4454
4455
4459 std::vector<Point<spacedim>> vertices;
4460
4464 std::vector<bool> vertices_used;
4465
4470 std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4472
4477
4478
4484
4495
4510 std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4512
4513
4533 std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4535
4536 // make a couple of classes friends
4537 template <int, int, int>
4538 friend class TriaAccessorBase;
4539 template <int, int, int>
4540 friend class TriaAccessor;
4541 friend class TriaAccessor<0, 1, spacedim>;
4542
4543 friend class CellAccessor<dim, spacedim>;
4544
4545 friend struct ::internal::TriaAccessorImplementation::Implementation;
4546
4547 friend struct ::internal::TriangulationImplementation::Implementation;
4548 friend struct ::internal::TriangulationImplementation::
4549 ImplementationMixedMesh;
4550
4551 friend class ::internal::TriangulationImplementation::TriaObjects;
4552
4553 // explicitly check for sensible template arguments, but not on windows
4554 // because MSVC creates bogus warnings during normal compilation
4555#ifndef DEAL_II_MSVC
4556 static_assert(dim <= spacedim,
4557 "The dimension <dim> of a Triangulation must be less than or "
4558 "equal to the space dimension <spacedim> in which it lives.");
4559#endif
4560};
4561
4562
4563#ifndef DOXYGEN
4564
4565
4566
4567namespace internal
4568{
4569 namespace TriangulationImplementation
4570 {
4571 template <class Archive>
4572 void
4573 NumberCache<1>::serialize(Archive &ar, const unsigned int)
4574 {
4575 ar &n_levels;
4576 ar &n_lines &n_lines_level;
4577 ar &n_active_lines &n_active_lines_level;
4578 }
4579
4580
4581 template <class Archive>
4582 void
4583 NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4584 {
4585 this->NumberCache<1>::serialize(ar, version);
4586
4587 ar &n_quads &n_quads_level;
4588 ar &n_active_quads &n_active_quads_level;
4589 }
4590
4591
4592 template <class Archive>
4593 void
4594 NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4595 {
4596 this->NumberCache<2>::serialize(ar, version);
4597
4598 ar &n_hexes &n_hexes_level;
4599 ar &n_active_hexes &n_active_hexes_level;
4600 }
4601
4602 } // namespace TriangulationImplementation
4603} // namespace internal
4604
4605
4606template <int dim, int spacedim>
4609 const unsigned int index) const
4610{
4611 AssertIndexRange(index, vertices_used.size());
4612 return vertices_used[index];
4613}
4614
4615
4616
4617template <int dim, int spacedim>
4619inline unsigned int Triangulation<dim, spacedim>::n_levels() const
4620{
4621 return number_cache.n_levels;
4622}
4623
4624template <int dim, int spacedim>
4626inline unsigned int Triangulation<dim, spacedim>::n_global_levels() const
4627{
4628 return number_cache.n_levels;
4629}
4630
4631
4632template <int dim, int spacedim>
4634inline unsigned int Triangulation<dim, spacedim>::n_vertices() const
4635{
4636 return vertices.size();
4637}
4638
4639
4640
4641template <int dim, int spacedim>
4643inline const std::vector<Point<spacedim>>
4645{
4646 return vertices;
4647}
4648
4649
4650template <int dim, int spacedim>
4652template <class Archive>
4653void Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4654{
4655 // as discussed in the documentation, do not store the signals as
4656 // well as boundary and manifold description but everything else
4657 ar &smooth_grid;
4658
4659 unsigned int n_levels = levels.size();
4660 ar &n_levels;
4661 for (const auto &level : levels)
4662 ar &level;
4663
4664 // boost dereferences a nullptr when serializing a nullptr
4665 // at least up to 1.65.1. This causes problems with clang-5.
4666 // Therefore, work around it.
4667 bool faces_is_nullptr = (faces.get() == nullptr);
4668 ar &faces_is_nullptr;
4669 if (!faces_is_nullptr)
4670 ar &faces;
4671
4672 ar &vertices;
4673 ar &vertices_used;
4674
4675 ar &anisotropic_refinement;
4676 ar &number_cache;
4677
4678 ar &check_for_distorted_cells;
4679
4680 if (dim == 1)
4681 {
4682 ar &vertex_to_boundary_id_map_1d;
4683 ar &vertex_to_manifold_id_map_1d;
4684 }
4685}
4686
4687
4688
4689template <int dim, int spacedim>
4691template <class Archive>
4692void Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4693{
4694 // clear previous content. this also calls the respective signal
4695 clear();
4696
4697 // as discussed in the documentation, do not store the signals as
4698 // well as boundary and manifold description but everything else
4699 ar &smooth_grid;
4700
4701 unsigned int size;
4702 ar &size;
4703 levels.resize(size);
4704 for (auto &level_ : levels)
4705 {
4706 std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4707 ar &level;
4708 level_ = std::move(level);
4709 }
4710
4711 // Workaround for nullptr, see in save().
4712 bool faces_is_nullptr = true;
4713 ar &faces_is_nullptr;
4714 if (!faces_is_nullptr)
4715 ar &faces;
4716
4717 ar &vertices;
4718 ar &vertices_used;
4719
4720 ar &anisotropic_refinement;
4721 ar &number_cache;
4722
4723 // the levels do not serialize the active_cell_indices because
4724 // they are easy enough to rebuild upon re-loading data. do
4725 // this here. don't forget to first resize the fields appropriately
4726 {
4727 for (auto &level : levels)
4728 {
4729 level->active_cell_indices.resize(level->refine_flags.size());
4730 level->global_active_cell_indices.resize(level->refine_flags.size());
4731 level->global_level_cell_indices.resize(level->refine_flags.size());
4732 }
4733 reset_cell_vertex_indices_cache();
4734 reset_active_cell_indices();
4735 reset_global_cell_indices();
4736 }
4737
4738 reset_policy();
4739
4740 bool my_check_for_distorted_cells;
4741 ar &my_check_for_distorted_cells;
4742
4743 Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4744 ExcMessage("The triangulation loaded into here must have the "
4745 "same setting with regard to reporting distorted "
4746 "cell as the one previously stored."));
4747
4748 if (dim == 1)
4749 {
4750 ar &vertex_to_boundary_id_map_1d;
4751 ar &vertex_to_manifold_id_map_1d;
4752 }
4753
4754 // trigger the create signal to indicate
4755 // that new content has been imported into
4756 // the triangulation
4757 signals.create();
4758}
4759
4760
4761
4762template <int dim, int spacedim>
4764inline unsigned int Triangulation<dim, spacedim>::
4766 const types::coarse_cell_id coarse_cell_id) const
4767{
4768 return coarse_cell_id;
4769}
4770
4771
4772
4773template <int dim, int spacedim>
4777 const unsigned int coarse_cell_index) const
4778{
4779 return coarse_cell_index;
4780}
4781
4782
4783
4784/* -------------- declaration of explicit specializations ------------- */
4785
4786template <>
4787unsigned int
4789template <>
4790unsigned int
4791Triangulation<1, 1>::n_quads(const unsigned int level) const;
4792template <>
4793unsigned int
4794Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4795template <>
4796unsigned int
4797Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4798template <>
4799unsigned int
4800Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4801template <>
4802unsigned int
4804template <>
4805unsigned int
4806Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4807template <>
4808unsigned int
4810template <>
4811unsigned int
4812Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4813template <>
4814unsigned int
4815Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4816template <>
4817unsigned int
4819template <>
4820unsigned int
4822template <>
4823unsigned int
4824Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4825template <>
4826unsigned int
4827Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4828
4829template <>
4830unsigned int
4832
4833
4834// -------------------------------------------------------------------
4835// -- Explicit specializations for codimension one grids
4836
4837
4838template <>
4839unsigned int
4841template <>
4842unsigned int
4843Triangulation<1, 2>::n_quads(const unsigned int level) const;
4844template <>
4845unsigned int
4846Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4847template <>
4848unsigned int
4849Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4850template <>
4851unsigned int
4852Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4853template <>
4854unsigned int
4855Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4856template <>
4857unsigned int
4859template <>
4860unsigned int
4862
4863// -------------------------------------------------------------------
4864// -- Explicit specializations for codimension two grids
4865
4866template <>
4867unsigned int
4869template <>
4870unsigned int
4871Triangulation<1, 3>::n_quads(const unsigned int level) const;
4872template <>
4873unsigned int
4874Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4875template <>
4876unsigned int
4877Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4878template <>
4879unsigned int
4880Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4881template <>
4882unsigned int
4883Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4884template <>
4885unsigned int
4887template <>
4888unsigned int
4890
4891template <>
4892bool
4894template <>
4895bool
4897template <>
4898bool
4900
4901
4902extern template class Triangulation<1, 1>;
4903extern template class Triangulation<1, 2>;
4904extern template class Triangulation<1, 3>;
4905extern template class Triangulation<2, 2>;
4906extern template class Triangulation<2, 3>;
4907extern template class Triangulation<3, 3>;
4908
4909#endif // DOXYGEN
4910
4912
4913// Include tria_accessor.h here, so that it is possible for an end
4914// user to use the iterators of Triangulation<dim> directly without
4915// the need to include tria_accessor.h separately. (Otherwise the
4916// iterators are an 'opaque' or 'incomplete' type.)
4918
4919#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:4088
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:4464
bool anisotropic_refinement
Definition tria.h:4476
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, unsigned char > > periodic_face_map
Definition tria.h:4073
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:4534
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:3906
const Triangulation< dim, spacedim > & get_triangulation() const
void reset_global_cell_indices()
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()
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:4459
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:4483
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:4511
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:4004
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:4453
::CellStatus CellStatus
Definition tria.h:2224
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:4089
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:4494
void load_user_flags_quad(const std::vector< bool > &v)
void load_user_flags_line(const std::vector< bool > &v)
void save_user_indices_hex(std::vector< unsigned int > &v) const
DistortedCellList execute_refinement()
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)
cell_iterator end() const
virtual bool has_hanging_nodes() const
std::vector< GridTools::PeriodicFacePair< cell_iterator > > periodic_face_pairs_level_0
Definition tria.h:4065
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
virtual void save(const std::string &filename) const
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:3856
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 save_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &filename) 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
virtual void load(const std::string &filename)
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
virtual void update_cell_relations()
Definition tria.h:3896
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
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:4090
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
Definition tria.h:4471
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition tria.h:3998
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition tria.h:4056
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition tria.h:2509
internal::CellAttachedDataSerializer< dim, spacedim > data_serializer
Definition tria.h:3908
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:4445
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 load_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const unsigned int local_num_cells, const std::string &filename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable)
void save_refine_flags(std::vector< bool > &v) const
std::vector< char > dest_data_variable
Definition tria.h:529
std::vector< char > src_data_fixed
Definition tria.h:519
std::vector< char > dest_data_fixed
Definition tria.h:520
std::vector< int > src_sizes_variable
Definition tria.h:526
std::vector< unsigned int > sizes_fixed_cumulative
Definition tria.h:513
typename std::pair< cell_iterator, CellStatus > cell_relation_t
Definition tria.h:387
std::vector< int > dest_sizes_variable
Definition tria.h:527
std::vector< char > src_data_variable
Definition tria.h:528
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
Point< 3 > vertices[4]
Point< 2 > first
Definition grid_out.cc:4623
unsigned int level
Definition grid_out.cc:4626
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:471
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516
static ::ExceptionBase & ExcMessage(std::string arg1)
typename IteratorSelector::hex_iterator hex_iterator
Definition tria.h:1685
typename IteratorSelector::active_quad_iterator active_quad_iterator
Definition tria.h:1676
typename IteratorSelector::active_hex_iterator active_hex_iterator
Definition tria.h:1696
typename IteratorSelector::quad_iterator quad_iterator
Definition tria.h:1661
typename IteratorSelector::line_iterator line_iterator
Definition tria.h:1637
typename IteratorSelector::active_line_iterator active_line_iterator
Definition tria.h:1652
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:2267
virtual ~DistortedCellList() noexcept override
boost::signals2::signal< void()> mesh_movement
Definition tria.h:2324
boost::signals2::signal< void()> post_distributed_load
Definition tria.h:2498
boost::signals2::signal< void()> pre_distributed_save
Definition tria.h:2476
boost::signals2::signal< unsigned int(const cell_iterator &, const ::CellStatus), CellWeightSum< unsigned int > > weight
Definition tria.h:2422
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition tria.h:2355
boost::signals2::signal< void()> pre_distributed_refinement
Definition tria.h:2435
boost::signals2::signal< void()> post_distributed_refinement
Definition tria.h:2454
boost::signals2::signal< void()> pre_refinement
Definition tria.h:2300
boost::signals2::signal< void()> any_change
Definition tria.h:2381
boost::signals2::signal< void()> create
Definition tria.h:2291
boost::signals2::signal< void()> clear
Definition tria.h:2370
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition tria.h:2345
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition tria.h:2335
boost::signals2::signal< void()> post_refinement
Definition tria.h:2307
boost::signals2::signal< void()> pre_partition
Definition tria.h:2315
boost::signals2::signal< void()> pre_distributed_repartition
Definition tria.h:2461
boost::signals2::signal< void()> pre_distributed_load
Definition tria.h:2491
boost::signals2::signal< void()> post_p4est_refinement
Definition tria.h:2444
boost::signals2::signal< void()> post_distributed_repartition
Definition tria.h:2468
boost::signals2::signal< void()> post_distributed_save
Definition tria.h:2483
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