Reference documentation for deal.II version GIT relicensing-1062-gc06da148b8 2024-07-15 19:20:02+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
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
1965 virtual std::vector<types::boundary_id>
1966 get_boundary_ids() const;
1967
1980 virtual std::vector<types::manifold_id>
1981 get_manifold_ids() const;
1982
2004 virtual void
2005 copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
2006
2055 virtual void
2056 create_triangulation(const std::vector<Point<spacedim>> &vertices,
2057 const std::vector<CellData<dim>> &cells,
2058 const SubCellData &subcelldata);
2059
2072 virtual void
2073 create_triangulation(
2074 const TriangulationDescription::Description<dim, spacedim>
2075 &construction_data);
2076
2085 void
2086 flip_all_direction_flags();
2087
2099 void
2100 set_all_refine_flags();
2101
2124 void
2125 refine_global(const unsigned int times = 1);
2126
2139 void
2140 coarsen_global(const unsigned int times = 1);
2141
2173 virtual void
2174 execute_coarsening_and_refinement();
2175
2205 virtual bool
2206 prepare_coarsening_and_refinement();
2207
2222 using CellStatus DEAL_II_DEPRECATED_EARLY = ::CellStatus;
2223
2228 static constexpr auto CELL_PERSIST DEAL_II_DEPRECATED_EARLY =
2230
2235 static constexpr auto CELL_REFINE DEAL_II_DEPRECATED_EARLY =
2237
2242 static constexpr auto CELL_COARSEN DEAL_II_DEPRECATED_EARLY =
2244
2249 static constexpr auto CELL_INVALID DEAL_II_DEPRECATED_EARLY =
2251
2252
2258 template <typename T>
2260 {
2261 using result_type = T;
2262
2263 template <typename InputIterator>
2264 T
2265 operator()(InputIterator first, InputIterator last) const
2266 {
2267 return std::accumulate(first, last, T());
2268 }
2269 };
2270
2280 struct Signals
2281 {
2289 boost::signals2::signal<void()> create;
2290
2298 boost::signals2::signal<void()> pre_refinement;
2299
2305 boost::signals2::signal<void()> post_refinement;
2306
2313 boost::signals2::signal<void()> pre_partition;
2314
2322 boost::signals2::signal<void()> mesh_movement;
2323
2331 boost::signals2::signal<void(
2332 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2334
2341 boost::signals2::signal<void(
2342 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2344
2351 boost::signals2::signal<void(
2352 const Triangulation<dim, spacedim> &destination_tria)>
2354
2368 boost::signals2::signal<void()> clear;
2369
2379 boost::signals2::signal<void()> any_change;
2380
2417 boost::signals2::signal<unsigned int(const cell_iterator &,
2418 const ::CellStatus),
2421
2433 boost::signals2::signal<void()> pre_distributed_refinement;
2434
2442 boost::signals2::signal<void()> post_p4est_refinement;
2443
2452 boost::signals2::signal<void()> post_distributed_refinement;
2453
2459 boost::signals2::signal<void()> pre_distributed_repartition;
2460
2466 boost::signals2::signal<void()> post_distributed_repartition;
2467
2474 boost::signals2::signal<void()> pre_distributed_save;
2475
2481 boost::signals2::signal<void()> post_distributed_save;
2482
2489 boost::signals2::signal<void()> pre_distributed_load;
2490
2496 boost::signals2::signal<void()> post_distributed_load;
2497 };
2498
2508
2520 void
2521 save_refine_flags(std::ostream &out) const;
2522
2526 void
2527 save_refine_flags(std::vector<bool> &v) const;
2528
2532 void
2533 load_refine_flags(std::istream &in);
2534
2538 void
2539 load_refine_flags(const std::vector<bool> &v);
2540
2544 void
2545 save_coarsen_flags(std::ostream &out) const;
2546
2550 void
2551 save_coarsen_flags(std::vector<bool> &v) const;
2552
2556 void
2557 load_coarsen_flags(std::istream &out);
2558
2562 void
2563 load_coarsen_flags(const std::vector<bool> &v);
2564
2569 bool
2571
2583 void
2585
2591 void
2592 save_user_flags(std::ostream &out) const;
2593
2599 void
2600 save_user_flags(std::vector<bool> &v) const;
2601
2606 void
2607 load_user_flags(std::istream &in);
2608
2613 void
2614 load_user_flags(const std::vector<bool> &v);
2615
2620 void
2622
2627 void
2628 save_user_flags_line(std::ostream &out) const;
2629
2635 void
2636 save_user_flags_line(std::vector<bool> &v) const;
2637
2642 void
2643 load_user_flags_line(std::istream &in);
2644
2649 void
2650 load_user_flags_line(const std::vector<bool> &v);
2651
2656 void
2658
2663 void
2664 save_user_flags_quad(std::ostream &out) const;
2665
2671 void
2672 save_user_flags_quad(std::vector<bool> &v) const;
2673
2678 void
2679 load_user_flags_quad(std::istream &in);
2680
2685 void
2686 load_user_flags_quad(const std::vector<bool> &v);
2687
2688
2693 void
2695
2700 void
2701 save_user_flags_hex(std::ostream &out) const;
2702
2708 void
2709 save_user_flags_hex(std::vector<bool> &v) const;
2710
2715 void
2716 load_user_flags_hex(std::istream &in);
2717
2722 void
2723 load_user_flags_hex(const std::vector<bool> &v);
2724
2730 void
2732
2738 void
2739 save_user_indices(std::vector<unsigned int> &v) const;
2740
2745 void
2746 load_user_indices(const std::vector<unsigned int> &v);
2747
2753 void
2754 save_user_pointers(std::vector<void *> &v) const;
2755
2760 void
2761 load_user_pointers(const std::vector<void *> &v);
2762
2768 void
2769 save_user_indices_line(std::vector<unsigned int> &v) const;
2770
2775 void
2776 load_user_indices_line(const std::vector<unsigned int> &v);
2777
2783 void
2784 save_user_indices_quad(std::vector<unsigned int> &v) const;
2785
2790 void
2791 load_user_indices_quad(const std::vector<unsigned int> &v);
2792
2798 void
2799 save_user_indices_hex(std::vector<unsigned int> &v) const;
2800
2805 void
2806 load_user_indices_hex(const std::vector<unsigned int> &v);
2812 void
2813 save_user_pointers_line(std::vector<void *> &v) const;
2814
2819 void
2820 load_user_pointers_line(const std::vector<void *> &v);
2821
2827 void
2828 save_user_pointers_quad(std::vector<void *> &v) const;
2829
2834 void
2835 load_user_pointers_quad(const std::vector<void *> &v);
2836
2842 void
2843 save_user_pointers_hex(std::vector<void *> &v) const;
2844
2849 void
2850 load_user_pointers_hex(const std::vector<void *> &v);
2851
2877 begin(const unsigned int level = 0) const;
2878
2910 begin_active(const unsigned int level = 0) const;
2911
2917 end() const;
2918
2938 end(const unsigned int level) const;
2939
2960 end_active(const unsigned int level) const;
2961
2962
2967 last() const;
2968
2974
2983 create_cell_iterator(const CellId &cell_id) const;
2984
2995 bool
2996 contains_cell(const CellId &cell_id) const;
3016
3054
3071 cell_iterators_on_level(const unsigned int level) const;
3072
3089 active_cell_iterators_on_level(const unsigned int level) const;
3090
3093 /*-------------------------------------------------------------------------*/
3094
3104 begin_face() const;
3105
3111
3117 end_face() const;
3118
3139
3142 /*-------------------------------------------------------------------------*/
3143
3155
3163
3170 end_vertex() const;
3171
3191 unsigned int
3192 n_lines() const;
3193
3197 unsigned int
3198 n_lines(const unsigned int level) const;
3199
3203 unsigned int
3205
3209 unsigned int
3210 n_active_lines(const unsigned int level) const;
3211
3215 unsigned int
3216 n_quads() const;
3217
3221 unsigned int
3222 n_quads(const unsigned int level) const;
3223
3227 unsigned int
3229
3233 unsigned int
3234 n_active_quads(const unsigned int level) const;
3235
3239 unsigned int
3240 n_hexs() const;
3241
3246 unsigned int
3247 n_hexs(const unsigned int level) const;
3248
3252 unsigned int
3254
3259 unsigned int
3260 n_active_hexs(const unsigned int level) const;
3261
3266 unsigned int
3267 n_cells() const;
3268
3273 unsigned int
3274 n_cells(const unsigned int level) const;
3275
3280 unsigned int
3282
3292
3293
3298 unsigned int
3299 n_active_cells(const unsigned int level) const;
3300
3305 virtual types::coarse_cell_id
3307
3313 unsigned int
3314 n_faces() const;
3315
3321 unsigned int
3323
3341 unsigned int
3342 n_levels() const;
3343
3350 virtual unsigned int
3352
3362 virtual bool
3364
3372 unsigned int
3373 n_vertices() const;
3374
3383 const std::vector<Point<spacedim>> &
3385
3390 unsigned int
3392
3396 bool
3397 vertex_used(const unsigned int index) const;
3398
3403 const std::vector<bool> &
3405
3417 unsigned int
3419
3426 virtual types::subdomain_id
3428
3440
3447
3448
3465 unsigned int
3467
3477 unsigned int
3478 n_raw_lines(const unsigned int level) const;
3479
3489 unsigned int
3491
3501 unsigned int
3502 n_raw_quads(const unsigned int level) const;
3503
3513 unsigned int
3514 n_raw_hexs(const unsigned int level) const;
3515
3525 unsigned int
3526 n_raw_cells(const unsigned int level) const;
3527
3539 unsigned int
3541
3552 virtual std::size_t
3554
3564 template <class Archive>
3565 void
3566 save(Archive &ar, const unsigned int version) const;
3567
3585 template <class Archive>
3586 void
3587 load(Archive &ar, const unsigned int version);
3588
3589
3596 virtual void
3597 save(const std::string &filename) const;
3598
3602 virtual void
3603 load(const std::string &filename);
3604
3605
3621 virtual void
3623 const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3624
3628 const std::map<
3629 std::pair<cell_iterator, unsigned int>,
3630 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>> &
3632
3637 const std::vector<ReferenceCell> &
3639
3644 bool
3646
3651 bool
3653
3659 bool
3661
3662#ifdef DOXYGEN
3668 template <class Archive>
3669 void
3670 serialize(Archive &archive, const unsigned int version);
3671#else
3672 // This macro defines the serialize() method that is compatible with
3673 // the templated save() and load() method that have been implemented.
3674 BOOST_SERIALIZATION_SPLIT_MEMBER()
3675#endif
3676
3681public:
3789 unsigned int
3791 const std::function<std::vector<char>(const cell_iterator &,
3792 const ::CellStatus)>
3793 &pack_callback,
3794 const bool returns_variable_size_data);
3795
3845 void
3847 const unsigned int handle,
3848 const std::function<
3849 void(const cell_iterator &,
3850 const ::CellStatus,
3851 const boost::iterator_range<std::vector<char>::const_iterator> &)>
3852 &unpack_callback);
3853
3855
3856protected:
3863 void
3864 save_attached_data(const unsigned int global_first_cell,
3865 const unsigned int global_num_cells,
3866 const std::string &filename) const;
3867
3875 void
3876 load_attached_data(const unsigned int global_first_cell,
3877 const unsigned int global_num_cells,
3878 const unsigned int local_num_cells,
3879 const std::string &filename,
3880 const unsigned int n_attached_deserialize_fixed,
3881 const unsigned int n_attached_deserialize_variable);
3882
3893 virtual void
3896
3902 std::vector<typename internal::CellAttachedDataSerializer<dim, spacedim>::
3903 cell_relation_t>
3905
3911public:
3922 DeclException2(ExcInvalidLevel,
3923 int,
3924 int,
3925 << "You are requesting information from refinement level "
3926 << arg1
3927 << " of a triangulation, but this triangulation only has "
3928 << arg2 << " refinement levels. The given level " << arg1
3929 << " must be *less* than " << arg2 << '.');
3937 ExcTriangulationNotEmpty,
3938 int,
3939 int,
3940 << "You are trying to perform an operation on a triangulation "
3941 << "that is only allowed if the triangulation is currently empty. "
3942 << "However, it currently stores " << arg1 << " vertices and has "
3943 << "cells on " << arg2 << " levels.");
3949 DeclException0(ExcGridReadError);
3954 DeclException0(ExcFacesHaveNoLevel);
3960 DeclException1(ExcEmptyLevel,
3961 int,
3962 << "You tried to do something on level " << arg1
3963 << ", but this level is empty.");
3964
3972 DeclException1(ExcBoundaryIdNotFound,
3974 << "The given boundary_id " << arg1
3975 << " is not defined in this Triangulation!");
3976
3983 ExcInconsistentCoarseningFlags,
3984 "A cell is flagged for coarsening, but either not all of its siblings "
3985 "are active or flagged for coarsening as well. Please clean up all "
3986 "coarsen flags on your triangulation via "
3987 "Triangulation::prepare_coarsening_and_refinement() beforehand!");
3988
3991protected:
3997
4002 std::vector<ReferenceCell> reference_cells;
4003
4017 static void
4018 write_bool_vector(const unsigned int magic_number1,
4019 const std::vector<bool> &v,
4020 const unsigned int magic_number2,
4021 std::ostream &out);
4022
4027 static void
4028 read_bool_vector(const unsigned int magic_number1,
4029 std::vector<bool> &v,
4030 const unsigned int magic_number2,
4031 std::istream &in);
4032
4037 void
4039
4043 virtual void
4045
4046
4047private:
4052 std::unique_ptr<
4055
4062 std::vector<GridTools::PeriodicFacePair<cell_iterator>>
4064
4069 std::map<std::pair<cell_iterator, unsigned int>,
4070 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>>
4072
4083 TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
4086 using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
4087 using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
4088 using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
4089
4100 begin_raw(const unsigned int level = 0) const;
4101
4107 end_raw(const unsigned int level) const;
4108
4123 begin_raw_line(const unsigned int level = 0) const;
4124
4143 begin_line(const unsigned int level = 0) const;
4144
4163 begin_active_line(const unsigned int level = 0) const;
4164
4170 end_line() const;
4171
4199 begin_raw_quad(const unsigned int level = 0) const;
4200
4219 begin_quad(const unsigned int level = 0) const;
4220
4239 begin_active_quad(const unsigned int level = 0) const;
4240
4246 end_quad() const;
4247
4274 begin_raw_hex(const unsigned int level = 0) const;
4275
4294 begin_hex(const unsigned int level = 0) const;
4295
4314 begin_active_hex(const unsigned int level = 0) const;
4315
4321 end_hex() const;
4322
4339 void
4341
4345 void
4347
4354 void
4356
4360 void
4362
4366 void
4368
4384
4391 void
4393
4398 void
4400
4416 virtual unsigned int
4418 const types::coarse_cell_id coarse_cell_id) const;
4419
4420
4433 virtual types::coarse_cell_id
4435 const unsigned int coarse_cell_index) const;
4436
4441 std::vector<
4442 std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
4444
4450 std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
4452
4453
4457 std::vector<Point<spacedim>> vertices;
4458
4462 std::vector<bool> vertices_used;
4463
4468 std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4470
4475
4476
4482
4493
4508 std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4510
4511
4531 std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4533
4534 // make a couple of classes friends
4535 template <int, int, int>
4536 friend class TriaAccessorBase;
4537 template <int, int, int>
4538 friend class TriaAccessor;
4539 friend class TriaAccessor<0, 1, spacedim>;
4540
4541 friend class CellAccessor<dim, spacedim>;
4542
4543 friend struct ::internal::TriaAccessorImplementation::Implementation;
4544
4545 friend struct ::internal::TriangulationImplementation::Implementation;
4546 friend struct ::internal::TriangulationImplementation::
4547 ImplementationMixedMesh;
4548
4549 friend class ::internal::TriangulationImplementation::TriaObjects;
4550
4551 // explicitly check for sensible template arguments, but not on windows
4552 // because MSVC creates bogus warnings during normal compilation
4553#ifndef DEAL_II_MSVC
4554 static_assert(dim <= spacedim,
4555 "The dimension <dim> of a Triangulation must be less than or "
4556 "equal to the space dimension <spacedim> in which it lives.");
4557#endif
4558};
4559
4560
4561#ifndef DOXYGEN
4562
4563
4564
4565namespace internal
4566{
4567 namespace TriangulationImplementation
4568 {
4569 template <class Archive>
4570 void
4571 NumberCache<1>::serialize(Archive &ar, const unsigned int)
4572 {
4573 ar &n_levels;
4574 ar &n_lines &n_lines_level;
4575 ar &n_active_lines &n_active_lines_level;
4576 }
4577
4578
4579 template <class Archive>
4580 void
4581 NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4582 {
4583 this->NumberCache<1>::serialize(ar, version);
4584
4585 ar &n_quads &n_quads_level;
4586 ar &n_active_quads &n_active_quads_level;
4587 }
4588
4589
4590 template <class Archive>
4591 void
4592 NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4593 {
4594 this->NumberCache<2>::serialize(ar, version);
4595
4596 ar &n_hexes &n_hexes_level;
4597 ar &n_active_hexes &n_active_hexes_level;
4598 }
4599
4600 } // namespace TriangulationImplementation
4601} // namespace internal
4602
4603
4604template <int dim, int spacedim>
4607 const unsigned int index) const
4608{
4609 AssertIndexRange(index, vertices_used.size());
4610 return vertices_used[index];
4611}
4612
4613
4614
4615template <int dim, int spacedim>
4617inline unsigned int Triangulation<dim, spacedim>::n_levels() const
4618{
4619 return number_cache.n_levels;
4620}
4621
4622template <int dim, int spacedim>
4624inline unsigned int Triangulation<dim, spacedim>::n_global_levels() const
4625{
4626 return number_cache.n_levels;
4627}
4628
4629
4630template <int dim, int spacedim>
4632inline unsigned int Triangulation<dim, spacedim>::n_vertices() const
4633{
4634 return vertices.size();
4635}
4636
4637
4638
4639template <int dim, int spacedim>
4641inline const std::vector<Point<spacedim>>
4643{
4644 return vertices;
4645}
4646
4647
4648template <int dim, int spacedim>
4650template <class Archive>
4651void Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4652{
4653 // as discussed in the documentation, do not store the signals as
4654 // well as boundary and manifold description but everything else
4655 ar &smooth_grid;
4656
4657 unsigned int n_levels = levels.size();
4658 ar &n_levels;
4659 for (const auto &level : levels)
4660 ar &level;
4661
4662 // boost dereferences a nullptr when serializing a nullptr
4663 // at least up to 1.65.1. This causes problems with clang-5.
4664 // Therefore, work around it.
4665 bool faces_is_nullptr = (faces.get() == nullptr);
4666 ar &faces_is_nullptr;
4667 if (!faces_is_nullptr)
4668 ar &faces;
4669
4670 ar &vertices;
4671 ar &vertices_used;
4672
4673 ar &anisotropic_refinement;
4674 ar &number_cache;
4675
4676 ar &check_for_distorted_cells;
4677
4678 if (dim == 1)
4679 {
4680 ar &vertex_to_boundary_id_map_1d;
4681 ar &vertex_to_manifold_id_map_1d;
4682 }
4683}
4684
4685
4686
4687template <int dim, int spacedim>
4689template <class Archive>
4690void Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4691{
4692 // clear previous content. this also calls the respective signal
4693 clear();
4694
4695 // as discussed in the documentation, do not store the signals as
4696 // well as boundary and manifold description but everything else
4697 ar &smooth_grid;
4698
4699 unsigned int size;
4700 ar &size;
4701 levels.resize(size);
4702 for (auto &level_ : levels)
4703 {
4704 std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4705 ar &level;
4706 level_ = std::move(level);
4707 }
4708
4709 // Workaround for nullptr, see in save().
4710 bool faces_is_nullptr = true;
4711 ar &faces_is_nullptr;
4712 if (!faces_is_nullptr)
4713 ar &faces;
4714
4715 ar &vertices;
4716 ar &vertices_used;
4717
4718 ar &anisotropic_refinement;
4719 ar &number_cache;
4720
4721 // the levels do not serialize the active_cell_indices because
4722 // they are easy enough to rebuild upon re-loading data. do
4723 // this here. don't forget to first resize the fields appropriately
4724 {
4725 for (const auto &level : levels)
4726 {
4727 level->active_cell_indices.resize(level->refine_flags.size());
4728 level->global_active_cell_indices.resize(level->refine_flags.size());
4729 level->global_level_cell_indices.resize(level->refine_flags.size());
4730 }
4731 reset_cell_vertex_indices_cache();
4732 reset_active_cell_indices();
4733 reset_global_cell_indices();
4734 }
4735
4736 reset_policy();
4737
4738 bool my_check_for_distorted_cells;
4739 ar &my_check_for_distorted_cells;
4740
4741 Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4742 ExcMessage("The triangulation loaded into here must have the "
4743 "same setting with regard to reporting distorted "
4744 "cell as the one previously stored."));
4745
4746 if (dim == 1)
4747 {
4748 ar &vertex_to_boundary_id_map_1d;
4749 ar &vertex_to_manifold_id_map_1d;
4750 }
4751
4752 // trigger the create signal to indicate
4753 // that new content has been imported into
4754 // the triangulation
4755 signals.create();
4756}
4757
4758
4759
4760template <int dim, int spacedim>
4762inline unsigned int Triangulation<dim, spacedim>::
4764 const types::coarse_cell_id coarse_cell_id) const
4765{
4766 return coarse_cell_id;
4767}
4768
4769
4770
4771template <int dim, int spacedim>
4775 const unsigned int coarse_cell_index) const
4776{
4777 return coarse_cell_index;
4778}
4779
4780
4781
4782/* -------------- declaration of explicit specializations ------------- */
4783
4784template <>
4785unsigned int
4787template <>
4788unsigned int
4789Triangulation<1, 1>::n_quads(const unsigned int level) const;
4790template <>
4791unsigned int
4792Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4793template <>
4794unsigned int
4795Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4796template <>
4797unsigned int
4798Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4799template <>
4800unsigned int
4802template <>
4803unsigned int
4804Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4805template <>
4806unsigned int
4808template <>
4809unsigned int
4810Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4811template <>
4812unsigned int
4813Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4814template <>
4815unsigned int
4817template <>
4818unsigned int
4820template <>
4821unsigned int
4822Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4823template <>
4824unsigned int
4825Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4826
4827template <>
4828unsigned int
4830
4831
4832// -------------------------------------------------------------------
4833// -- Explicit specializations for codimension one grids
4834
4835
4836template <>
4837unsigned int
4839template <>
4840unsigned int
4841Triangulation<1, 2>::n_quads(const unsigned int level) const;
4842template <>
4843unsigned int
4844Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4845template <>
4846unsigned int
4847Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4848template <>
4849unsigned int
4850Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4851template <>
4852unsigned int
4853Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4854template <>
4855unsigned int
4857template <>
4858unsigned int
4860
4861// -------------------------------------------------------------------
4862// -- Explicit specializations for codimension two grids
4863
4864template <>
4865unsigned int
4867template <>
4868unsigned int
4869Triangulation<1, 3>::n_quads(const unsigned int level) const;
4870template <>
4871unsigned int
4872Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4873template <>
4874unsigned int
4875Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4876template <>
4877unsigned int
4878Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4879template <>
4880unsigned int
4881Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4882template <>
4883unsigned int
4885template <>
4886unsigned int
4888
4889template <>
4890bool
4892template <>
4893bool
4895template <>
4896bool
4898
4899
4900extern template class Triangulation<1, 1>;
4901extern template class Triangulation<1, 2>;
4902extern template class Triangulation<1, 3>;
4903extern template class Triangulation<2, 2>;
4904extern template class Triangulation<2, 3>;
4905extern template class Triangulation<3, 3>;
4906
4907#endif // DOXYGEN
4908
4910
4911// Include tria_accessor.h here, so that it is possible for an end
4912// user to use the iterators of Triangulation<dim> directly without
4913// the need to include tria_accessor.h separately. (Otherwise the
4914// iterators are an 'opaque' or 'incomplete' type.)
4916
4917#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:4086
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:4462
bool anisotropic_refinement
Definition tria.h:4474
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, unsigned char > > periodic_face_map
Definition tria.h:4071
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:4532
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:3904
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:4457
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:4481
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:4509
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:4002
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:4451
::CellStatus CellStatus
Definition tria.h:2222
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:4087
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:4492
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:4063
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:3854
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:3894
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:4088
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
Definition tria.h:4469
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition tria.h:3996
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition tria.h:4054
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition tria.h:2507
internal::CellAttachedDataSerializer< dim, spacedim > data_serializer
Definition tria.h:3906
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:4443
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:2265
virtual ~DistortedCellList() noexcept override
boost::signals2::signal< void()> mesh_movement
Definition tria.h:2322
boost::signals2::signal< void()> post_distributed_load
Definition tria.h:2496
boost::signals2::signal< void()> pre_distributed_save
Definition tria.h:2474
boost::signals2::signal< unsigned int(const cell_iterator &, const ::CellStatus), CellWeightSum< unsigned int > > weight
Definition tria.h:2420
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition tria.h:2353
boost::signals2::signal< void()> pre_distributed_refinement
Definition tria.h:2433
boost::signals2::signal< void()> post_distributed_refinement
Definition tria.h:2452
boost::signals2::signal< void()> pre_refinement
Definition tria.h:2298
boost::signals2::signal< void()> any_change
Definition tria.h:2379
boost::signals2::signal< void()> create
Definition tria.h:2289
boost::signals2::signal< void()> clear
Definition tria.h:2368
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition tria.h:2343
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition tria.h:2333
boost::signals2::signal< void()> post_refinement
Definition tria.h:2305
boost::signals2::signal< void()> pre_partition
Definition tria.h:2313
boost::signals2::signal< void()> pre_distributed_repartition
Definition tria.h:2459
boost::signals2::signal< void()> pre_distributed_load
Definition tria.h:2489
boost::signals2::signal< void()> post_p4est_refinement
Definition tria.h:2442
boost::signals2::signal< void()> post_distributed_repartition
Definition tria.h:2466
boost::signals2::signal< void()> post_distributed_save
Definition tria.h:2481
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