deal.II version GIT relicensing-2850-g1646a780ac 2025-03-17 15:10:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
tria.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_tria_h
16#define dealii_tria_h
17
18
19#include <deal.II/base/config.h>
20
26#include <deal.II/base/point.h>
27
33
34#include <boost/range/iterator_range_core.hpp>
35#include <boost/serialization/map.hpp>
36#include <boost/serialization/split_member.hpp>
37#include <boost/serialization/unique_ptr.hpp>
38#include <boost/serialization/vector.hpp>
39#include <boost/signals2.hpp>
40
41#include <bitset>
42#include <functional>
43#include <list>
44#include <map>
45#include <memory>
46#include <numeric>
47#include <vector>
48
49
51
52#ifdef signals
53# error \
54 "The name 'signals' is already defined. You are most likely using the QT library \
55and using the 'signals' keyword. You can either #include the Qt headers (or any conflicting headers) \
56*after* the deal.II headers or you can define the 'QT_NO_KEYWORDS' macro and use the 'Q_SIGNALS' macro."
57#endif
58
59// Forward declarations
60#ifndef DOXYGEN
61template <int dim, int spacedim>
62class Manifold;
63
64template <int dim>
65struct CellData;
66
67struct SubCellData;
68
70{
71 template <int, int>
72 struct Description;
73}
74
75namespace GridTools
76{
77 template <typename CellIterator>
78 struct PeriodicFacePair;
79}
80
81template <int, int, int>
82class TriaAccessor;
83template <int spacedim>
84class TriaAccessor<0, 1, spacedim>;
85template <int, int, int>
87
88namespace internal
89{
90 namespace TriangulationImplementation
91 {
92 class TriaFaces;
93
94 class TriaObjects;
95
96 template <int, int>
97 class Policy;
98
104 struct Implementation;
105 struct ImplementationMixedMesh;
106 } // namespace TriangulationImplementation
107
108 namespace TriaAccessorImplementation
109 {
110 struct Implementation;
111 }
112} // namespace internal
113#endif
114
115
116/*------------------------------------------------------------------------*/
117
118
119namespace internal
120{
125 namespace TriangulationImplementation
126 {
138 template <int dim>
140 {};
141
153 template <>
154 struct NumberCache<1>
155 {
159 unsigned int n_levels;
160
164 unsigned int n_lines;
165
169 std::vector<unsigned int> n_lines_level;
170
174 unsigned int n_active_lines;
175
179 std::vector<unsigned int> n_active_lines_level;
180
184 std::shared_ptr<const Utilities::MPI::Partitioner>
186
190 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
192
196 NumberCache();
197
202 std::size_t
203 memory_consumption() const;
204
210 template <class Archive>
211 void
212 serialize(Archive &ar, const unsigned int version);
213 };
214
215
228 template <>
229 struct NumberCache<2> : public NumberCache<1>
230 {
234 unsigned int n_quads;
235
239 std::vector<unsigned int> n_quads_level;
240
244 unsigned int n_active_quads;
245
249 std::vector<unsigned int> n_active_quads_level;
250
254 NumberCache();
255
260 std::size_t
261 memory_consumption() const;
262
268 template <class Archive>
269 void
270 serialize(Archive &ar, const unsigned int version);
271 };
272
273
287 template <>
288 struct NumberCache<3> : public NumberCache<2>
289 {
293 unsigned int n_hexes;
294
298 std::vector<unsigned int> n_hexes_level;
299
303 unsigned int n_active_hexes;
304
308 std::vector<unsigned int> n_active_hexes_level;
309
313 NumberCache();
314
319 std::size_t
320 memory_consumption() const;
321
327 template <class Archive>
328 void
329 serialize(Archive &ar, const unsigned int version);
330 };
331 } // namespace TriangulationImplementation
332
333
337 template <int dim, int spacedim = dim>
340 {
342
348
354
356 std::function<std::vector<char>(cell_iterator, CellStatus)>;
357
362 std::vector<pack_callback_t> pack_callbacks_fixed;
363 std::vector<pack_callback_t> pack_callbacks_variable;
364 };
365
376 template <int dim, int spacedim = dim>
379 {
380 public:
382
387 static inline constexpr unsigned int version_number = 5;
388
394 using cell_relation_t = typename std::pair<cell_iterator, CellStatus>;
395
397
406 void
407 pack_data(
408 const std::vector<cell_relation_t> &cell_relations,
409 const std::vector<
411 &pack_callbacks_fixed,
412 const std::vector<
414 &pack_callbacks_variable,
415 const MPI_Comm &mpi_communicator);
416
424 void
425 unpack_cell_status(std::vector<cell_relation_t> &cell_relations) const;
426
439 void
440 unpack_data(
441 const std::vector<cell_relation_t> &cell_relations,
442 const unsigned int handle,
443 const std::function<
444 void(const cell_iterator &,
445 const CellStatus &,
446 const boost::iterator_range<std::vector<char>::const_iterator> &)>
447 &unpack_callback) const;
448
463 void
464 save(const unsigned int global_first_cell,
465 const unsigned int global_num_cells,
466 const std::string &file_basename,
467 const MPI_Comm &mpi_communicator) const;
468
487 void
488 load(const unsigned int global_first_cell,
489 const unsigned int global_num_cells,
490 const unsigned int local_num_cells,
491 const std::string &file_basename,
492 const unsigned int n_attached_deserialize_fixed,
493 const unsigned int n_attached_deserialize_variable,
494 const MPI_Comm &mpi_communicator);
495
502 void
503 clear();
504
509
520 std::vector<unsigned int> sizes_fixed_cumulative;
521
526 std::vector<char> src_data_fixed;
527 std::vector<char> dest_data_fixed;
528
533 std::vector<int> src_sizes_variable;
534 std::vector<int> dest_sizes_variable;
535 std::vector<char> src_data_variable;
536 std::vector<char> dest_data_variable;
537 };
538} // namespace internal
539
540
541/*------------------------------------------------------------------------*/
542
543
1321template <int dim, int spacedim = dim>
1324{
1325private:
1332
1333public:
1339 {
1344 none = 0x0,
1388 limit_level_difference_at_vertices = 0x1,
1409 eliminate_unrefined_islands = 0x2,
1425 patch_level_1 = 0x4,
1449 coarsest_level_1 = 0x8,
1474 allow_anisotropic_smoothing = 0x10,
1507 eliminate_refined_inner_islands = 0x100,
1512 eliminate_refined_boundary_islands = 0x200,
1518 do_not_produce_unrefined_islands = 0x400,
1519
1524 smoothing_on_refinement =
1525 (limit_level_difference_at_vertices | eliminate_unrefined_islands),
1530 smoothing_on_coarsening =
1531 (eliminate_refined_inner_islands | eliminate_refined_boundary_islands |
1532 do_not_produce_unrefined_islands),
1533
1539 maximum_smoothing = 0xffff ^ allow_anisotropic_smoothing
1541
1558
1564
1582
1596 using face_iterator = TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1597
1610 TriaActiveIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1611
1621
1636
1644 using line_iterator = typename IteratorSelector::line_iterator;
1645
1659 using active_line_iterator = typename IteratorSelector::active_line_iterator;
1660
1668 using quad_iterator = typename IteratorSelector::quad_iterator;
1669
1683 using active_quad_iterator = typename IteratorSelector::active_quad_iterator;
1684
1692 using hex_iterator = typename IteratorSelector::hex_iterator;
1693
1703 using active_hex_iterator = typename IteratorSelector::active_hex_iterator;
1704
1724 {
1731 virtual ~DistortedCellList() noexcept override;
1732
1737 std::list<typename Triangulation<dim, spacedim>::cell_iterator>
1738 distorted_cells;
1739 };
1740
1744 static constexpr unsigned int dimension = dim;
1745
1749 static constexpr unsigned int space_dimension = spacedim;
1750
1765 Triangulation(const MeshSmoothing smooth_grid = none,
1766 const bool check_for_distorted_cells = false);
1767
1783 Triangulation(const Triangulation<dim, spacedim> &) = delete;
1784
1791 Triangulation(Triangulation<dim, spacedim> &&tria) noexcept;
1792
1797 operator=(Triangulation<dim, spacedim> &&tria) noexcept;
1798
1802 virtual ~Triangulation() override;
1803
1810 virtual void
1811 clear();
1812
1817 virtual MPI_Comm
1818 get_mpi_communicator() const;
1819
1826 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
1827 "Access the MPI communicator with get_mpi_communicator() instead.")
1828 MPI_Comm
1829 get_communicator() const;
1830
1836 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1837 global_active_cell_index_partitioner() const;
1838
1844 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1845 global_level_cell_index_partitioner(const unsigned int level) const;
1846
1851 virtual void
1852 set_mesh_smoothing(const MeshSmoothing mesh_smoothing);
1853
1857 virtual const MeshSmoothing &
1858 get_mesh_smoothing() const;
1859
1882 void
1883 set_manifold(const types::manifold_id number,
1884 const Manifold<dim, spacedim> &manifold_object);
1885
1902 void
1903 reset_manifold(const types::manifold_id manifold_number);
1904
1921 void
1922 reset_all_manifolds();
1923
1932 void
1933 set_all_manifold_ids(const types::manifold_id number);
1934
1943 void
1944 set_all_manifold_ids_on_boundary(const types::manifold_id number);
1945
1955 void
1956 set_all_manifold_ids_on_boundary(const types::boundary_id b_id,
1957 const types::manifold_id number);
1958
1972 const Manifold<dim, spacedim> &
1973 get_manifold(const types::manifold_id number) const;
1974
1985 virtual std::vector<types::boundary_id>
1986 get_boundary_ids() const;
1987
2000 virtual std::vector<types::manifold_id>
2001 get_manifold_ids() const;
2002
2024 virtual void
2025 copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
2026
2075 virtual void
2076 create_triangulation(const std::vector<Point<spacedim>> &vertices,
2077 const std::vector<CellData<dim>> &cells,
2078 const SubCellData &subcelldata);
2079
2092 virtual void
2093 create_triangulation(
2094 const TriangulationDescription::Description<dim, spacedim>
2095 &construction_data);
2096
2105 void
2106 flip_all_direction_flags();
2107
2119 void
2120 set_all_refine_flags();
2121
2144 void
2145 refine_global(const unsigned int times = 1);
2146
2159 void
2160 coarsen_global(const unsigned int times = 1);
2161
2193 virtual void
2194 execute_coarsening_and_refinement();
2195
2225 virtual bool
2226 prepare_coarsening_and_refinement();
2227
2243
2248 static constexpr auto CELL_PERSIST DEAL_II_DEPRECATED =
2250
2255 static constexpr auto CELL_REFINE DEAL_II_DEPRECATED =
2257
2262 static constexpr auto CELL_COARSEN DEAL_II_DEPRECATED =
2264
2269 static constexpr auto CELL_INVALID DEAL_II_DEPRECATED =
2271
2272
2278 template <typename T>
2280 {
2281 using result_type = T;
2282
2283 template <typename InputIterator>
2284 T
2285 operator()(InputIterator first, InputIterator last) const
2286 {
2287 return std::accumulate(first, last, T());
2288 }
2289 };
2290
2300 struct Signals
2301 {
2309 boost::signals2::signal<void()> create;
2310
2318 boost::signals2::signal<void()> pre_refinement;
2319
2325 boost::signals2::signal<void()> post_refinement;
2326
2333 boost::signals2::signal<void()> pre_partition;
2334
2342 boost::signals2::signal<void()> mesh_movement;
2343
2351 boost::signals2::signal<void(
2352 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2354
2361 boost::signals2::signal<void(
2362 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2364
2371 boost::signals2::signal<void(
2372 const Triangulation<dim, spacedim> &destination_tria)>
2374
2388 boost::signals2::signal<void()> clear;
2389
2399 boost::signals2::signal<void()> any_change;
2400
2437 boost::signals2::signal<unsigned int(const cell_iterator &,
2438 const ::CellStatus),
2441
2453 boost::signals2::signal<void()> pre_distributed_refinement;
2454
2462 boost::signals2::signal<void()> post_p4est_refinement;
2463
2472 boost::signals2::signal<void()> post_distributed_refinement;
2473
2479 boost::signals2::signal<void()> pre_distributed_repartition;
2480
2486 boost::signals2::signal<void()> post_distributed_repartition;
2487
2494 boost::signals2::signal<void()> pre_distributed_save;
2495
2501 boost::signals2::signal<void()> post_distributed_save;
2502
2509 boost::signals2::signal<void()> pre_distributed_load;
2510
2516 boost::signals2::signal<void()> post_distributed_load;
2517 };
2518
2528
2540 void
2541 save_refine_flags(std::ostream &out) const;
2542
2546 void
2547 save_refine_flags(std::vector<bool> &v) const;
2548
2552 void
2553 load_refine_flags(std::istream &in);
2554
2558 void
2559 load_refine_flags(const std::vector<bool> &v);
2560
2564 void
2565 save_coarsen_flags(std::ostream &out) const;
2566
2570 void
2571 save_coarsen_flags(std::vector<bool> &v) const;
2572
2576 void
2577 load_coarsen_flags(std::istream &out);
2578
2582 void
2583 load_coarsen_flags(const std::vector<bool> &v);
2584
2589 bool
2591
2603 void
2605
2611 void
2612 save_user_flags(std::ostream &out) const;
2613
2619 void
2620 save_user_flags(std::vector<bool> &v) const;
2621
2626 void
2627 load_user_flags(std::istream &in);
2628
2633 void
2634 load_user_flags(const std::vector<bool> &v);
2635
2640 void
2642
2647 void
2648 save_user_flags_line(std::ostream &out) const;
2649
2655 void
2656 save_user_flags_line(std::vector<bool> &v) const;
2657
2662 void
2663 load_user_flags_line(std::istream &in);
2664
2669 void
2670 load_user_flags_line(const std::vector<bool> &v);
2671
2676 void
2678
2683 void
2684 save_user_flags_quad(std::ostream &out) const;
2685
2691 void
2692 save_user_flags_quad(std::vector<bool> &v) const;
2693
2698 void
2699 load_user_flags_quad(std::istream &in);
2700
2705 void
2706 load_user_flags_quad(const std::vector<bool> &v);
2707
2708
2713 void
2715
2720 void
2721 save_user_flags_hex(std::ostream &out) const;
2722
2728 void
2729 save_user_flags_hex(std::vector<bool> &v) const;
2730
2735 void
2736 load_user_flags_hex(std::istream &in);
2737
2742 void
2743 load_user_flags_hex(const std::vector<bool> &v);
2744
2750 void
2752
2758 void
2759 save_user_indices(std::vector<unsigned int> &v) const;
2760
2765 void
2766 load_user_indices(const std::vector<unsigned int> &v);
2767
2773 void
2774 save_user_pointers(std::vector<void *> &v) const;
2775
2780 void
2781 load_user_pointers(const std::vector<void *> &v);
2782
2788 void
2789 save_user_indices_line(std::vector<unsigned int> &v) const;
2790
2795 void
2796 load_user_indices_line(const std::vector<unsigned int> &v);
2797
2803 void
2804 save_user_indices_quad(std::vector<unsigned int> &v) const;
2805
2810 void
2811 load_user_indices_quad(const std::vector<unsigned int> &v);
2812
2818 void
2819 save_user_indices_hex(std::vector<unsigned int> &v) const;
2820
2825 void
2826 load_user_indices_hex(const std::vector<unsigned int> &v);
2832 void
2833 save_user_pointers_line(std::vector<void *> &v) const;
2834
2839 void
2840 load_user_pointers_line(const std::vector<void *> &v);
2841
2847 void
2848 save_user_pointers_quad(std::vector<void *> &v) const;
2849
2854 void
2855 load_user_pointers_quad(const std::vector<void *> &v);
2856
2862 void
2863 save_user_pointers_hex(std::vector<void *> &v) const;
2864
2869 void
2870 load_user_pointers_hex(const std::vector<void *> &v);
2871
2897 begin(const unsigned int level = 0) const;
2898
2930 begin_active(const unsigned int level = 0) const;
2931
2937 end() const;
2938
2958 end(const unsigned int level) const;
2959
2980 end_active(const unsigned int level) const;
2981
2982
2987 last() const;
2988
2994
3003 create_cell_iterator(const CellId &cell_id) const;
3004
3015 bool
3016 contains_cell(const CellId &cell_id) const;
3036
3074
3091 cell_iterators_on_level(const unsigned int level) const;
3092
3109 active_cell_iterators_on_level(const unsigned int level) const;
3110
3113 /*-------------------------------------------------------------------------*/
3114
3124 begin_face() const;
3125
3131
3137 end_face() const;
3138
3159
3162 /*-------------------------------------------------------------------------*/
3163
3175
3183
3190 end_vertex() const;
3191
3211 unsigned int
3212 n_lines() const;
3213
3217 unsigned int
3218 n_lines(const unsigned int level) const;
3219
3223 unsigned int
3225
3229 unsigned int
3230 n_active_lines(const unsigned int level) const;
3231
3235 unsigned int
3236 n_quads() const;
3237
3241 unsigned int
3242 n_quads(const unsigned int level) const;
3243
3247 unsigned int
3249
3253 unsigned int
3254 n_active_quads(const unsigned int level) const;
3255
3259 unsigned int
3260 n_hexs() const;
3261
3266 unsigned int
3267 n_hexs(const unsigned int level) const;
3268
3272 unsigned int
3274
3279 unsigned int
3280 n_active_hexs(const unsigned int level) const;
3281
3286 unsigned int
3287 n_cells() const;
3288
3293 unsigned int
3294 n_cells(const unsigned int level) const;
3295
3300 unsigned int
3302
3312
3313
3318 unsigned int
3319 n_active_cells(const unsigned int level) const;
3320
3325 virtual types::coarse_cell_id
3327
3333 unsigned int
3334 n_faces() const;
3335
3341 unsigned int
3343
3361 unsigned int
3362 n_levels() const;
3363
3370 virtual unsigned int
3372
3382 virtual bool
3384
3392 unsigned int
3393 n_vertices() const;
3394
3403 const std::vector<Point<spacedim>> &
3405
3410 unsigned int
3412
3416 bool
3417 vertex_used(const unsigned int index) const;
3418
3423 const std::vector<bool> &
3425
3437 unsigned int
3439
3446 virtual types::subdomain_id
3448
3460
3467
3468
3485 unsigned int
3487
3497 unsigned int
3498 n_raw_lines(const unsigned int level) const;
3499
3509 unsigned int
3511
3521 unsigned int
3522 n_raw_quads(const unsigned int level) const;
3523
3533 unsigned int
3534 n_raw_hexs(const unsigned int level) const;
3535
3545 unsigned int
3546 n_raw_cells(const unsigned int level) const;
3547
3559 unsigned int
3561
3572 virtual std::size_t
3574
3584 template <class Archive>
3585 void
3586 save(Archive &ar, const unsigned int version) const;
3587
3605 template <class Archive>
3606 void
3607 load(Archive &ar, const unsigned int version);
3608
3609
3618 virtual void
3619 save(const std::string &file_basename) const;
3620
3624 virtual void
3625 load(const std::string &file_basename);
3626
3627
3643 virtual void
3645 const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3646
3650 const std::map<std::pair<cell_iterator, unsigned int>,
3651 std::pair<std::pair<cell_iterator, unsigned int>,
3654
3659 const std::vector<ReferenceCell> &
3661
3666 bool
3668
3673 bool
3675
3681 bool
3683
3684#ifdef DOXYGEN
3692 template <class Archive>
3693 void
3694 serialize(Archive &archive, const unsigned int version);
3695#else
3696 // This macro defines the serialize() method that is compatible with
3697 // the templated save() and load() method that have been implemented.
3698 BOOST_SERIALIZATION_SPLIT_MEMBER()
3699#endif
3700
3705public:
3813 unsigned int
3815 const std::function<std::vector<char>(const cell_iterator &,
3816 const ::CellStatus)>
3817 &pack_callback,
3818 const bool returns_variable_size_data);
3819
3869 void
3871 const unsigned int handle,
3872 const std::function<
3873 void(const cell_iterator &,
3874 const ::CellStatus,
3875 const boost::iterator_range<std::vector<char>::const_iterator> &)>
3876 &unpack_callback);
3877
3879
3880protected:
3888 void
3889 save_attached_data(const unsigned int global_first_cell,
3890 const unsigned int global_num_cells,
3891 const std::string &file_basename) const;
3892
3901 void
3902 load_attached_data(const unsigned int global_first_cell,
3903 const unsigned int global_num_cells,
3904 const unsigned int local_num_cells,
3905 const std::string &file_basename,
3906 const unsigned int n_attached_deserialize_fixed,
3907 const unsigned int n_attached_deserialize_variable);
3908
3918 void
3920
3926 void
3928
3929
3934 void
3936
3942 std::vector<typename internal::CellAttachedDataSerializer<dim, spacedim>::
3943 cell_relation_t>
3945
3951public:
3962 DeclException2(ExcInvalidLevel,
3963 int,
3964 int,
3965 << "You are requesting information from refinement level "
3966 << arg1
3967 << " of a triangulation, but this triangulation only has "
3968 << arg2 << " refinement levels. The given level " << arg1
3969 << " must be *less* than " << arg2 << '.');
3977 ExcTriangulationNotEmpty,
3978 int,
3979 int,
3980 << "You are trying to perform an operation on a triangulation "
3981 << "that is only allowed if the triangulation is currently empty. "
3982 << "However, it currently stores " << arg1 << " vertices and has "
3983 << "cells on " << arg2 << " levels.");
3989 DeclException0(ExcGridReadError);
3994 DeclException0(ExcFacesHaveNoLevel);
4000 DeclException1(ExcEmptyLevel,
4001 int,
4002 << "You tried to do something on level " << arg1
4003 << ", but this level is empty.");
4004
4012 DeclException1(ExcBoundaryIdNotFound,
4014 << "The given boundary_id " << arg1
4015 << " is not defined in this Triangulation!");
4016
4023 ExcInconsistentCoarseningFlags,
4024 "A cell is flagged for coarsening, but either not all of its siblings "
4025 "are active or flagged for coarsening as well. Please clean up all "
4026 "coarsen flags on your triangulation via "
4027 "Triangulation::prepare_coarsening_and_refinement() beforehand!");
4028
4031protected:
4037
4042 std::vector<ReferenceCell> reference_cells;
4043
4057 static void
4058 write_bool_vector(const unsigned int magic_number1,
4059 const std::vector<bool> &v,
4060 const unsigned int magic_number2,
4061 std::ostream &out);
4062
4067 static void
4068 read_bool_vector(const unsigned int magic_number1,
4069 std::vector<bool> &v,
4070 const unsigned int magic_number2,
4071 std::istream &in);
4072
4077 void
4079
4083 virtual void
4085
4086
4087private:
4092 std::unique_ptr<
4095
4102 std::vector<GridTools::PeriodicFacePair<cell_iterator>>
4104
4109 std::map<std::pair<cell_iterator, unsigned int>,
4110 std::pair<std::pair<cell_iterator, unsigned int>,
4113
4124 TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
4127 using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
4128 using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
4129 using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
4130
4141 begin_raw(const unsigned int level = 0) const;
4142
4148 end_raw(const unsigned int level) const;
4149
4164 begin_raw_line(const unsigned int level = 0) const;
4165
4184 begin_line(const unsigned int level = 0) const;
4185
4204 begin_active_line(const unsigned int level = 0) const;
4205
4211 end_line() const;
4212
4240 begin_raw_quad(const unsigned int level = 0) const;
4241
4260 begin_quad(const unsigned int level = 0) const;
4261
4280 begin_active_quad(const unsigned int level = 0) const;
4281
4287 end_quad() const;
4288
4315 begin_raw_hex(const unsigned int level = 0) const;
4316
4335 begin_hex(const unsigned int level = 0) const;
4336
4355 begin_active_hex(const unsigned int level = 0) const;
4356
4362 end_hex() const;
4363
4380 void
4382
4386 void
4388
4395 void
4397
4401 void
4403
4407 void
4409
4425
4432 void
4434
4439 void
4441
4457 virtual unsigned int
4459 const types::coarse_cell_id coarse_cell_id) const;
4460
4461
4474 virtual types::coarse_cell_id
4476 const unsigned int coarse_cell_index) const;
4477
4482 std::vector<
4483 std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
4485
4491 std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
4493
4494
4498 std::vector<Point<spacedim>> vertices;
4499
4503 std::vector<bool> vertices_used;
4504
4508 std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4510
4515
4516
4522
4533
4548 std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4550
4551
4571 std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4573
4574 // make a couple of classes friends
4575 template <int, int, int>
4576 friend class TriaAccessorBase;
4577 template <int, int, int>
4578 friend class TriaAccessor;
4579 friend class TriaAccessor<0, 1, spacedim>;
4580
4581 friend class CellAccessor<dim, spacedim>;
4582
4583 friend struct ::internal::TriaAccessorImplementation::Implementation;
4584
4585 friend struct ::internal::TriangulationImplementation::Implementation;
4586 friend struct ::internal::TriangulationImplementation::
4587 ImplementationMixedMesh;
4588
4589 friend class ::internal::TriangulationImplementation::TriaObjects;
4590
4591 // explicitly check for sensible template arguments, but not on windows
4592 // because MSVC creates bogus warnings during normal compilation
4593#ifndef DEAL_II_MSVC
4594 static_assert(dim <= spacedim,
4595 "The dimension <dim> of a Triangulation must be less than or "
4596 "equal to the space dimension <spacedim> in which it lives.");
4597#endif
4598};
4599
4600
4601#ifndef DOXYGEN
4602
4603
4604
4605namespace internal
4606{
4607 namespace TriangulationImplementation
4608 {
4609 template <class Archive>
4610 void
4611 NumberCache<1>::serialize(Archive &ar, const unsigned int)
4612 {
4613 ar &n_levels;
4614 ar &n_lines &n_lines_level;
4615 ar &n_active_lines &n_active_lines_level;
4616 }
4617
4618
4619 template <class Archive>
4620 void
4621 NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4622 {
4623 this->NumberCache<1>::serialize(ar, version);
4624
4625 ar &n_quads &n_quads_level;
4626 ar &n_active_quads &n_active_quads_level;
4627 }
4628
4629
4630 template <class Archive>
4631 void
4632 NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4633 {
4634 this->NumberCache<2>::serialize(ar, version);
4635
4636 ar &n_hexes &n_hexes_level;
4637 ar &n_active_hexes &n_active_hexes_level;
4638 }
4639
4640 } // namespace TriangulationImplementation
4641} // namespace internal
4642
4643
4644template <int dim, int spacedim>
4647 const unsigned int index) const
4648{
4649 AssertIndexRange(index, vertices_used.size());
4650 return vertices_used[index];
4651}
4652
4653
4654
4655template <int dim, int spacedim>
4657inline unsigned int Triangulation<dim, spacedim>::n_levels() const
4658{
4659 return number_cache.n_levels;
4660}
4661
4662template <int dim, int spacedim>
4664inline unsigned int Triangulation<dim, spacedim>::n_global_levels() const
4665{
4666 return number_cache.n_levels;
4667}
4668
4669
4670template <int dim, int spacedim>
4672inline unsigned int Triangulation<dim, spacedim>::n_vertices() const
4673{
4674 return vertices.size();
4675}
4676
4677
4678
4679template <int dim, int spacedim>
4681inline const std::vector<Point<spacedim>>
4683{
4684 return vertices;
4685}
4686
4687
4688template <int dim, int spacedim>
4690template <class Archive>
4691void Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4692{
4693 // as discussed in the documentation, do not store the signals as
4694 // well as boundary and manifold description but everything else
4695 ar &smooth_grid;
4696
4697 unsigned int n_levels = levels.size();
4698 ar &n_levels;
4699 for (const auto &level : levels)
4700 ar &level;
4701
4702 // boost dereferences a nullptr when serializing a nullptr
4703 // at least up to 1.65.1. This causes problems with clang-5.
4704 // Therefore, work around it.
4705 bool faces_is_nullptr = (faces.get() == nullptr);
4706 ar &faces_is_nullptr;
4707 if (!faces_is_nullptr)
4708 ar &faces;
4709
4710 ar &vertices;
4711 ar &vertices_used;
4712
4713 ar &anisotropic_refinement;
4714 ar &number_cache;
4715
4716 ar &check_for_distorted_cells;
4717
4718 if (dim == 1)
4719 {
4720 ar &vertex_to_boundary_id_map_1d;
4721 ar &vertex_to_manifold_id_map_1d;
4722 }
4723}
4724
4725
4726
4727template <int dim, int spacedim>
4729template <class Archive>
4730void Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4731{
4732 // clear previous content. this also calls the respective signal
4733 clear();
4734
4735 // as discussed in the documentation, do not store the signals as
4736 // well as boundary and manifold description but everything else
4737 ar &smooth_grid;
4738
4739 unsigned int size;
4740 ar &size;
4741 levels.resize(size);
4742 for (auto &level_ : levels)
4743 {
4744 std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4745 ar &level;
4746 level_ = std::move(level);
4747 }
4748
4749 // Workaround for nullptr, see in save().
4750 bool faces_is_nullptr = true;
4751 ar &faces_is_nullptr;
4752 if (!faces_is_nullptr)
4753 ar &faces;
4754
4755 ar &vertices;
4756 ar &vertices_used;
4757
4758 ar &anisotropic_refinement;
4759 ar &number_cache;
4760
4761 // the levels do not serialize the active_cell_indices because
4762 // they are easy enough to rebuild upon re-loading data. do
4763 // this here. don't forget to first resize the fields appropriately
4764 {
4765 for (const auto &level : levels)
4766 {
4767 level->active_cell_indices.resize(level->refine_flags.size());
4768 level->global_active_cell_indices.resize(level->refine_flags.size());
4769 level->global_level_cell_indices.resize(level->refine_flags.size());
4770 }
4771 reset_cell_vertex_indices_cache();
4772 reset_active_cell_indices();
4773 reset_global_cell_indices();
4774 }
4775
4776 reset_policy();
4777
4778 bool my_check_for_distorted_cells;
4779 ar &my_check_for_distorted_cells;
4780
4781 Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4782 ExcMessage("The triangulation loaded into here must have the "
4783 "same setting with regard to reporting distorted "
4784 "cell as the one previously stored."));
4785
4786 if (dim == 1)
4787 {
4788 ar &vertex_to_boundary_id_map_1d;
4789 ar &vertex_to_manifold_id_map_1d;
4790 }
4791
4792 // trigger the create signal to indicate
4793 // that new content has been imported into
4794 // the triangulation
4795 signals.create();
4796}
4797
4798
4799
4800template <int dim, int spacedim>
4802inline unsigned int Triangulation<dim, spacedim>::
4804 const types::coarse_cell_id coarse_cell_id) const
4805{
4806 return coarse_cell_id;
4807}
4808
4809
4810
4811template <int dim, int spacedim>
4815 const unsigned int coarse_cell_index) const
4816{
4817 return coarse_cell_index;
4818}
4819
4820
4821
4822/* -------------- declaration of explicit specializations ------------- */
4823
4824template <>
4825unsigned int
4827template <>
4828unsigned int
4829Triangulation<1, 1>::n_quads(const unsigned int level) const;
4830template <>
4831unsigned int
4832Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4833template <>
4834unsigned int
4835Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4836template <>
4837unsigned int
4838Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4839template <>
4840unsigned int
4842template <>
4843unsigned int
4844Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4845template <>
4846unsigned int
4848template <>
4849unsigned int
4850Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4851template <>
4852unsigned int
4853Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4854template <>
4855unsigned int
4857template <>
4858unsigned int
4860template <>
4861unsigned int
4862Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4863template <>
4864unsigned int
4865Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4866
4867template <>
4868unsigned int
4870
4871
4872// -------------------------------------------------------------------
4873// -- Explicit specializations for codimension one grids
4874
4875
4876template <>
4877unsigned int
4879template <>
4880unsigned int
4881Triangulation<1, 2>::n_quads(const unsigned int level) const;
4882template <>
4883unsigned int
4884Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4885template <>
4886unsigned int
4887Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4888template <>
4889unsigned int
4890Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4891template <>
4892unsigned int
4893Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4894template <>
4895unsigned int
4897template <>
4898unsigned int
4900
4901// -------------------------------------------------------------------
4902// -- Explicit specializations for codimension two grids
4903
4904template <>
4905unsigned int
4907template <>
4908unsigned int
4909Triangulation<1, 3>::n_quads(const unsigned int level) const;
4910template <>
4911unsigned int
4912Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4913template <>
4914unsigned int
4915Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4916template <>
4917unsigned int
4918Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4919template <>
4920unsigned int
4921Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4922template <>
4923unsigned int
4925template <>
4926unsigned int
4928
4929template <>
4930bool
4932template <>
4933bool
4935template <>
4936bool
4938
4939
4940extern template class Triangulation<1, 1>;
4941extern template class Triangulation<1, 2>;
4942extern template class Triangulation<1, 3>;
4943extern template class Triangulation<2, 2>;
4944extern template class Triangulation<2, 3>;
4945extern template class Triangulation<3, 3>;
4946
4947#endif // DOXYGEN
4948
4950
4951#endif
CellStatus
Definition cell_status.h:31
@ cell_will_be_refined
@ children_will_be_coarsened
Definition point.h:113
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:4127
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:4503
bool anisotropic_refinement
Definition tria.h:4514
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:4572
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
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:3944
const Triangulation< dim, spacedim > & get_triangulation() const
void reset_global_cell_indices()
virtual void save(const std::string &file_basename) const
face_iterator end_face() const
void reset_active_cell_indices()
cell_iterator create_cell_iterator(const CellId &cell_id) const
cell_iterator begin(const unsigned int level=0) const
void fix_coarsen_flags()
virtual void load(const std::string &file_basename)
unsigned int n_active_cells(const unsigned int level) const
unsigned int n_cells(const unsigned int level) const
void save_user_pointers_line(std::vector< void * > &v) const
void save_user_flags_quad(std::vector< bool > &v) const
void load_refine_flags(std::istream &in)
void save_user_indices_line(std::vector< unsigned int > &v) const
raw_cell_iterator begin_raw(const unsigned int level=0) const
unsigned int n_lines() const
void load_coarsen_flags(const std::vector< bool > &v)
unsigned int n_raw_lines() const
virtual std::size_t memory_consumption() const
std::vector< Point< spacedim > > vertices
Definition tria.h:4498
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:4521
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:4549
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:4042
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:4492
::CellStatus CellStatus
Definition tria.h:2242
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:4128
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:4532
void load_user_flags_quad(const std::vector< bool > &v)
void save_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &file_basename) const
void load_user_flags_line(const std::vector< bool > &v)
void save_user_indices_hex(std::vector< unsigned int > &v) const
DistortedCellList execute_refinement()
void update_cell_relations()
active_line_iterator begin_active_line(const unsigned int level=0) const
void save_user_indices_quad(std::vector< unsigned int > &v) const
void load_user_pointers_hex(const std::vector< void * > &v)
void pack_data_serial()
cell_iterator end() const
virtual bool has_hanging_nodes() const
std::vector< GridTools::PeriodicFacePair< cell_iterator > > periodic_face_pairs_level_0
Definition tria.h:4103
unsigned int n_raw_cells(const unsigned int level) const
cell_iterator end(const unsigned int level) const
bool contains_cell(const CellId &cell_id) const
void load_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const unsigned int local_num_cells, const std::string &file_basename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable)
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const
unsigned int n_raw_quads(const unsigned int level) const
internal::CellAttachedData< dim, spacedim > cell_attached_data
Definition tria.h:3878
void load_coarsen_flags(std::istream &out)
quad_iterator end_quad() const
line_iterator begin_line(const unsigned int level=0) const
unsigned int max_adjacent_cells() const
vertex_iterator begin_vertex() const
void clear_user_flags()
unsigned int n_hexs() const
vertex_iterator end_vertex() const
bool vertex_used(const unsigned int index) const
void load_user_pointers_line(const std::vector< void * > &v)
void save_user_flags(std::vector< bool > &v) const
hex_iterator end_hex() const
hex_iterator begin_hex(const unsigned int level=0) const
virtual unsigned int n_global_levels() const
active_cell_iterator end_active(const unsigned int level) const
bool is_mixed_mesh() const
cell_iterator last() const
unsigned int n_active_quads() const
void save_coarsen_flags(std::vector< bool > &v) const
void load_user_indices_hex(const std::vector< unsigned int > &v)
unsigned int n_raw_quads() const
void save_user_pointers(std::vector< void * > &v) const
face_iterator begin_face() const
unsigned int n_cells() const
virtual bool prepare_coarsening_and_refinement()
unsigned int n_active_quads(const unsigned int level) const
void unpack_data_serial()
const std::vector< bool > & get_used_vertices() const
void load_user_flags(const std::vector< bool > &v)
typename IteratorSelector::raw_hex_iterator raw_hex_iterator
Definition tria.h:4129
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
Definition tria.h:4509
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition tria.h:4036
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition tria.h:4094
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition tria.h:2527
internal::CellAttachedDataSerializer< dim, spacedim > data_serializer
Definition tria.h:3946
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, types::geometric_orientation > > periodic_face_map
Definition tria.h:4112
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:4484
unsigned int n_raw_hexs(const unsigned int level) const
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, types::geometric_orientation > > & get_periodic_face_map() const
unsigned int n_lines(const unsigned int level) const
unsigned int n_active_hexs() const
virtual unsigned int coarse_cell_id_to_coarse_cell_index(const types::coarse_cell_id coarse_cell_id) const
unsigned int n_quads(const unsigned int level) const
void load_user_flags(std::istream &in)
void reset_policy()
void save_coarsen_flags(std::ostream &out) const
active_face_iterator begin_active_face() const
void clear_user_flags_line()
raw_line_iterator begin_raw_line(const unsigned int level=0) const
static void write_bool_vector(const unsigned int magic_number1, const std::vector< bool > &v, const unsigned int magic_number2, std::ostream &out)
active_cell_iterator begin_active(const unsigned int level=0) const
void execute_coarsening()
void save_refine_flags(std::vector< bool > &v) const
std::vector< char > dest_data_variable
Definition tria.h:536
std::vector< char > src_data_fixed
Definition tria.h:526
std::vector< char > dest_data_fixed
Definition tria.h:527
std::vector< int > src_sizes_variable
Definition tria.h:533
std::vector< unsigned int > sizes_fixed_cumulative
Definition tria.h:520
typename std::pair< cell_iterator, CellStatus > cell_relation_t
Definition tria.h:394
std::vector< int > dest_sizes_variable
Definition tria.h:534
std::vector< char > src_data_variable
Definition tria.h:535
#define DEAL_II_DEPRECATED
Definition config.h:283
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:245
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
Point< 2 > first
Definition grid_out.cc:4632
unsigned int level
Definition grid_out.cc:4635
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)
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcMessage(std::string arg1)
typename IteratorSelector::hex_iterator hex_iterator
Definition tria.h:1692
typename IteratorSelector::active_quad_iterator active_quad_iterator
Definition tria.h:1683
typename IteratorSelector::active_hex_iterator active_hex_iterator
Definition tria.h:1703
typename IteratorSelector::quad_iterator quad_iterator
Definition tria.h:1668
typename IteratorSelector::line_iterator line_iterator
Definition tria.h:1644
typename IteratorSelector::active_line_iterator active_line_iterator
Definition tria.h:1659
std::size_t size
Definition mpi.cc:745
STL namespace.
Definition types.h:32
global_cell_index coarse_cell_id
Definition types.h:139
unsigned char geometric_orientation
Definition types.h:40
T operator()(InputIterator first, InputIterator last) const
Definition tria.h:2285
virtual ~DistortedCellList() noexcept override
boost::signals2::signal< void()> mesh_movement
Definition tria.h:2342
boost::signals2::signal< void()> post_distributed_load
Definition tria.h:2516
boost::signals2::signal< void()> pre_distributed_save
Definition tria.h:2494
boost::signals2::signal< unsigned int(const cell_iterator &, const ::CellStatus), CellWeightSum< unsigned int > > weight
Definition tria.h:2440
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition tria.h:2373
boost::signals2::signal< void()> pre_distributed_refinement
Definition tria.h:2453
boost::signals2::signal< void()> post_distributed_refinement
Definition tria.h:2472
boost::signals2::signal< void()> pre_refinement
Definition tria.h:2318
boost::signals2::signal< void()> any_change
Definition tria.h:2399
boost::signals2::signal< void()> create
Definition tria.h:2309
boost::signals2::signal< void()> clear
Definition tria.h:2388
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition tria.h:2363
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition tria.h:2353
boost::signals2::signal< void()> post_refinement
Definition tria.h:2325
boost::signals2::signal< void()> pre_partition
Definition tria.h:2333
boost::signals2::signal< void()> pre_distributed_repartition
Definition tria.h:2479
boost::signals2::signal< void()> pre_distributed_load
Definition tria.h:2509
boost::signals2::signal< void()> post_p4est_refinement
Definition tria.h:2462
boost::signals2::signal< void()> post_distributed_repartition
Definition tria.h:2486
boost::signals2::signal< void()> post_distributed_save
Definition tria.h:2501
std::vector< pack_callback_t > pack_callbacks_variable
Definition tria.h:363
unsigned int n_attached_deserialize
Definition tria.h:353
std::vector< pack_callback_t > pack_callbacks_fixed
Definition tria.h:362
unsigned int n_attached_data_sets
Definition tria.h:347
std::function< std::vector< char >(cell_iterator, CellStatus)> pack_callback_t
Definition tria.h:356
std::shared_ptr< const Utilities::MPI::Partitioner > active_cell_index_partitioner
Definition tria.h:185
void serialize(Archive &ar, const unsigned int version)
std::vector< unsigned int > n_active_lines_level
Definition tria.h:179
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > level_cell_index_partitioners
Definition tria.h:191
std::vector< unsigned int > n_active_quads_level
Definition tria.h:249
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:308