Reference documentation for deal.II version 9.6.0
\(\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
386 static inline constexpr unsigned int version_number = 5;
387
393 using cell_relation_t = typename std::pair<cell_iterator, CellStatus>;
394
396
405 void
406 pack_data(
407 const std::vector<cell_relation_t> &cell_relations,
408 const std::vector<
410 &pack_callbacks_fixed,
411 const std::vector<
413 &pack_callbacks_variable,
414 const MPI_Comm &mpi_communicator);
415
423 void
424 unpack_cell_status(std::vector<cell_relation_t> &cell_relations) const;
425
438 void
439 unpack_data(
440 const std::vector<cell_relation_t> &cell_relations,
441 const unsigned int handle,
442 const std::function<
443 void(const cell_iterator &,
444 const CellStatus &,
445 const boost::iterator_range<std::vector<char>::const_iterator> &)>
446 &unpack_callback) const;
447
462 void
463 save(const unsigned int global_first_cell,
464 const unsigned int global_num_cells,
465 const std::string &file_basename,
466 const MPI_Comm &mpi_communicator) const;
467
486 void
487 load(const unsigned int global_first_cell,
488 const unsigned int global_num_cells,
489 const unsigned int local_num_cells,
490 const std::string &file_basename,
491 const unsigned int n_attached_deserialize_fixed,
492 const unsigned int n_attached_deserialize_variable,
493 const MPI_Comm &mpi_communicator);
494
501 void
502 clear();
503
508
519 std::vector<unsigned int> sizes_fixed_cumulative;
520
525 std::vector<char> src_data_fixed;
526 std::vector<char> dest_data_fixed;
527
532 std::vector<int> src_sizes_variable;
533 std::vector<int> dest_sizes_variable;
534 std::vector<char> src_data_variable;
535 std::vector<char> dest_data_variable;
536 };
537} // namespace internal
538
539
540/*------------------------------------------------------------------------*/
541
542
1320template <int dim, int spacedim = dim>
1323{
1324private:
1331
1332public:
1338 {
1343 none = 0x0,
1387 limit_level_difference_at_vertices = 0x1,
1408 eliminate_unrefined_islands = 0x2,
1424 patch_level_1 = 0x4,
1448 coarsest_level_1 = 0x8,
1473 allow_anisotropic_smoothing = 0x10,
1506 eliminate_refined_inner_islands = 0x100,
1511 eliminate_refined_boundary_islands = 0x200,
1517 do_not_produce_unrefined_islands = 0x400,
1518
1523 smoothing_on_refinement =
1524 (limit_level_difference_at_vertices | eliminate_unrefined_islands),
1529 smoothing_on_coarsening =
1530 (eliminate_refined_inner_islands | eliminate_refined_boundary_islands |
1531 do_not_produce_unrefined_islands),
1532
1538 maximum_smoothing = 0xffff ^ allow_anisotropic_smoothing
1540
1557
1563
1581
1595 using face_iterator = TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1596
1609 TriaActiveIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1610
1620
1635
1643 using line_iterator = typename IteratorSelector::line_iterator;
1644
1658 using active_line_iterator = typename IteratorSelector::active_line_iterator;
1659
1667 using quad_iterator = typename IteratorSelector::quad_iterator;
1668
1682 using active_quad_iterator = typename IteratorSelector::active_quad_iterator;
1683
1691 using hex_iterator = typename IteratorSelector::hex_iterator;
1692
1702 using active_hex_iterator = typename IteratorSelector::active_hex_iterator;
1703
1723 {
1730 virtual ~DistortedCellList() noexcept override;
1731
1736 std::list<typename Triangulation<dim, spacedim>::cell_iterator>
1737 distorted_cells;
1738 };
1739
1743 static constexpr unsigned int dimension = dim;
1744
1748 static constexpr unsigned int space_dimension = spacedim;
1749
1764 Triangulation(const MeshSmoothing smooth_grid = none,
1765 const bool check_for_distorted_cells = false);
1766
1782 Triangulation(const Triangulation<dim, spacedim> &) = delete;
1783
1790 Triangulation(Triangulation<dim, spacedim> &&tria) noexcept;
1791
1796 operator=(Triangulation<dim, spacedim> &&tria) noexcept;
1797
1801 virtual ~Triangulation() override;
1802
1809 virtual void
1810 clear();
1811
1816 virtual MPI_Comm
1817 get_communicator() const;
1818
1824 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1825 global_active_cell_index_partitioner() const;
1826
1832 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1833 global_level_cell_index_partitioner(const unsigned int level) const;
1834
1839 virtual void
1840 set_mesh_smoothing(const MeshSmoothing mesh_smoothing);
1841
1845 virtual const MeshSmoothing &
1846 get_mesh_smoothing() const;
1847
1870 void
1871 set_manifold(const types::manifold_id number,
1872 const Manifold<dim, spacedim> &manifold_object);
1873
1890 void
1891 reset_manifold(const types::manifold_id manifold_number);
1892
1909 void
1910 reset_all_manifolds();
1911
1920 void
1921 set_all_manifold_ids(const types::manifold_id number);
1922
1931 void
1932 set_all_manifold_ids_on_boundary(const types::manifold_id number);
1933
1943 void
1944 set_all_manifold_ids_on_boundary(const types::boundary_id b_id,
1945 const types::manifold_id number);
1946
1958 const Manifold<dim, spacedim> &
1959 get_manifold(const types::manifold_id number) const;
1960
1971 virtual std::vector<types::boundary_id>
1972 get_boundary_ids() const;
1973
1986 virtual std::vector<types::manifold_id>
1987 get_manifold_ids() const;
1988
2010 virtual void
2011 copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
2012
2061 virtual void
2062 create_triangulation(const std::vector<Point<spacedim>> &vertices,
2063 const std::vector<CellData<dim>> &cells,
2064 const SubCellData &subcelldata);
2065
2078 virtual void
2079 create_triangulation(
2080 const TriangulationDescription::Description<dim, spacedim>
2081 &construction_data);
2082
2091 void
2092 flip_all_direction_flags();
2093
2105 void
2106 set_all_refine_flags();
2107
2130 void
2131 refine_global(const unsigned int times = 1);
2132
2145 void
2146 coarsen_global(const unsigned int times = 1);
2147
2179 virtual void
2180 execute_coarsening_and_refinement();
2181
2211 virtual bool
2212 prepare_coarsening_and_refinement();
2213
2228 using CellStatus DEAL_II_DEPRECATED_EARLY = ::CellStatus;
2229
2234 static constexpr auto CELL_PERSIST DEAL_II_DEPRECATED_EARLY =
2236
2241 static constexpr auto CELL_REFINE DEAL_II_DEPRECATED_EARLY =
2243
2248 static constexpr auto CELL_COARSEN DEAL_II_DEPRECATED_EARLY =
2250
2255 static constexpr auto CELL_INVALID DEAL_II_DEPRECATED_EARLY =
2257
2258
2264 template <typename T>
2266 {
2267 using result_type = T;
2268
2269 template <typename InputIterator>
2270 T
2271 operator()(InputIterator first, InputIterator last) const
2272 {
2273 return std::accumulate(first, last, T());
2274 }
2275 };
2276
2286 struct Signals
2287 {
2295 boost::signals2::signal<void()> create;
2296
2304 boost::signals2::signal<void()> pre_refinement;
2305
2311 boost::signals2::signal<void()> post_refinement;
2312
2319 boost::signals2::signal<void()> pre_partition;
2320
2328 boost::signals2::signal<void()> mesh_movement;
2329
2337 boost::signals2::signal<void(
2338 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2340
2347 boost::signals2::signal<void(
2348 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2350
2357 boost::signals2::signal<void(
2358 const Triangulation<dim, spacedim> &destination_tria)>
2360
2374 boost::signals2::signal<void()> clear;
2375
2385 boost::signals2::signal<void()> any_change;
2386
2423 boost::signals2::signal<unsigned int(const cell_iterator &,
2424 const ::CellStatus),
2427
2439 boost::signals2::signal<void()> pre_distributed_refinement;
2440
2448 boost::signals2::signal<void()> post_p4est_refinement;
2449
2458 boost::signals2::signal<void()> post_distributed_refinement;
2459
2465 boost::signals2::signal<void()> pre_distributed_repartition;
2466
2472 boost::signals2::signal<void()> post_distributed_repartition;
2473
2480 boost::signals2::signal<void()> pre_distributed_save;
2481
2487 boost::signals2::signal<void()> post_distributed_save;
2488
2495 boost::signals2::signal<void()> pre_distributed_load;
2496
2502 boost::signals2::signal<void()> post_distributed_load;
2503 };
2504
2514
2526 void
2527 save_refine_flags(std::ostream &out) const;
2528
2532 void
2533 save_refine_flags(std::vector<bool> &v) const;
2534
2538 void
2539 load_refine_flags(std::istream &in);
2540
2544 void
2545 load_refine_flags(const std::vector<bool> &v);
2546
2550 void
2551 save_coarsen_flags(std::ostream &out) const;
2552
2556 void
2557 save_coarsen_flags(std::vector<bool> &v) const;
2558
2562 void
2563 load_coarsen_flags(std::istream &out);
2564
2568 void
2569 load_coarsen_flags(const std::vector<bool> &v);
2570
2575 bool
2577
2589 void
2591
2597 void
2598 save_user_flags(std::ostream &out) const;
2599
2605 void
2606 save_user_flags(std::vector<bool> &v) const;
2607
2612 void
2613 load_user_flags(std::istream &in);
2614
2619 void
2620 load_user_flags(const std::vector<bool> &v);
2621
2626 void
2628
2633 void
2634 save_user_flags_line(std::ostream &out) const;
2635
2641 void
2642 save_user_flags_line(std::vector<bool> &v) const;
2643
2648 void
2649 load_user_flags_line(std::istream &in);
2650
2655 void
2656 load_user_flags_line(const std::vector<bool> &v);
2657
2662 void
2664
2669 void
2670 save_user_flags_quad(std::ostream &out) const;
2671
2677 void
2678 save_user_flags_quad(std::vector<bool> &v) const;
2679
2684 void
2685 load_user_flags_quad(std::istream &in);
2686
2691 void
2692 load_user_flags_quad(const std::vector<bool> &v);
2693
2694
2699 void
2701
2706 void
2707 save_user_flags_hex(std::ostream &out) const;
2708
2714 void
2715 save_user_flags_hex(std::vector<bool> &v) const;
2716
2721 void
2722 load_user_flags_hex(std::istream &in);
2723
2728 void
2729 load_user_flags_hex(const std::vector<bool> &v);
2730
2736 void
2738
2744 void
2745 save_user_indices(std::vector<unsigned int> &v) const;
2746
2751 void
2752 load_user_indices(const std::vector<unsigned int> &v);
2753
2759 void
2760 save_user_pointers(std::vector<void *> &v) const;
2761
2766 void
2767 load_user_pointers(const std::vector<void *> &v);
2768
2774 void
2775 save_user_indices_line(std::vector<unsigned int> &v) const;
2776
2781 void
2782 load_user_indices_line(const std::vector<unsigned int> &v);
2783
2789 void
2790 save_user_indices_quad(std::vector<unsigned int> &v) const;
2791
2796 void
2797 load_user_indices_quad(const std::vector<unsigned int> &v);
2798
2804 void
2805 save_user_indices_hex(std::vector<unsigned int> &v) const;
2806
2811 void
2812 load_user_indices_hex(const std::vector<unsigned int> &v);
2818 void
2819 save_user_pointers_line(std::vector<void *> &v) const;
2820
2825 void
2826 load_user_pointers_line(const std::vector<void *> &v);
2827
2833 void
2834 save_user_pointers_quad(std::vector<void *> &v) const;
2835
2840 void
2841 load_user_pointers_quad(const std::vector<void *> &v);
2842
2848 void
2849 save_user_pointers_hex(std::vector<void *> &v) const;
2850
2855 void
2856 load_user_pointers_hex(const std::vector<void *> &v);
2857
2883 begin(const unsigned int level = 0) const;
2884
2916 begin_active(const unsigned int level = 0) const;
2917
2923 end() const;
2924
2944 end(const unsigned int level) const;
2945
2966 end_active(const unsigned int level) const;
2967
2968
2973 last() const;
2974
2980
2989 create_cell_iterator(const CellId &cell_id) const;
2990
3001 bool
3002 contains_cell(const CellId &cell_id) const;
3022
3060
3077 cell_iterators_on_level(const unsigned int level) const;
3078
3095 active_cell_iterators_on_level(const unsigned int level) const;
3096
3099 /*-------------------------------------------------------------------------*/
3100
3110 begin_face() const;
3111
3117
3123 end_face() const;
3124
3145
3148 /*-------------------------------------------------------------------------*/
3149
3161
3169
3176 end_vertex() const;
3177
3197 unsigned int
3198 n_lines() const;
3199
3203 unsigned int
3204 n_lines(const unsigned int level) const;
3205
3209 unsigned int
3211
3215 unsigned int
3216 n_active_lines(const unsigned int level) const;
3217
3221 unsigned int
3222 n_quads() const;
3223
3227 unsigned int
3228 n_quads(const unsigned int level) const;
3229
3233 unsigned int
3235
3239 unsigned int
3240 n_active_quads(const unsigned int level) const;
3241
3245 unsigned int
3246 n_hexs() const;
3247
3252 unsigned int
3253 n_hexs(const unsigned int level) const;
3254
3258 unsigned int
3260
3265 unsigned int
3266 n_active_hexs(const unsigned int level) const;
3267
3272 unsigned int
3273 n_cells() const;
3274
3279 unsigned int
3280 n_cells(const unsigned int level) const;
3281
3286 unsigned int
3288
3298
3299
3304 unsigned int
3305 n_active_cells(const unsigned int level) const;
3306
3311 virtual types::coarse_cell_id
3313
3319 unsigned int
3320 n_faces() const;
3321
3327 unsigned int
3329
3347 unsigned int
3348 n_levels() const;
3349
3356 virtual unsigned int
3358
3368 virtual bool
3370
3378 unsigned int
3379 n_vertices() const;
3380
3389 const std::vector<Point<spacedim>> &
3391
3396 unsigned int
3398
3402 bool
3403 vertex_used(const unsigned int index) const;
3404
3409 const std::vector<bool> &
3411
3423 unsigned int
3425
3432 virtual types::subdomain_id
3434
3446
3453
3454
3471 unsigned int
3473
3483 unsigned int
3484 n_raw_lines(const unsigned int level) const;
3485
3495 unsigned int
3497
3507 unsigned int
3508 n_raw_quads(const unsigned int level) const;
3509
3519 unsigned int
3520 n_raw_hexs(const unsigned int level) const;
3521
3531 unsigned int
3532 n_raw_cells(const unsigned int level) const;
3533
3545 unsigned int
3547
3558 virtual std::size_t
3560
3570 template <class Archive>
3571 void
3572 save(Archive &ar, const unsigned int version) const;
3573
3591 template <class Archive>
3592 void
3593 load(Archive &ar, const unsigned int version);
3594
3595
3604 virtual void
3605 save(const std::string &file_basename) const;
3606
3610 virtual void
3611 load(const std::string &file_basename);
3612
3613
3629 virtual void
3631 const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3632
3636 const std::map<
3637 std::pair<cell_iterator, unsigned int>,
3638 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>> &
3640
3645 const std::vector<ReferenceCell> &
3647
3652 bool
3654
3659 bool
3661
3667 bool
3669
3670#ifdef DOXYGEN
3678 template <class Archive>
3679 void
3680 serialize(Archive &archive, const unsigned int version);
3681#else
3682 // This macro defines the serialize() method that is compatible with
3683 // the templated save() and load() method that have been implemented.
3684 BOOST_SERIALIZATION_SPLIT_MEMBER()
3685#endif
3686
3691public:
3799 unsigned int
3801 const std::function<std::vector<char>(const cell_iterator &,
3802 const ::CellStatus)>
3803 &pack_callback,
3804 const bool returns_variable_size_data);
3805
3855 void
3857 const unsigned int handle,
3858 const std::function<
3859 void(const cell_iterator &,
3860 const ::CellStatus,
3861 const boost::iterator_range<std::vector<char>::const_iterator> &)>
3862 &unpack_callback);
3863
3865
3866protected:
3874 void
3875 save_attached_data(const unsigned int global_first_cell,
3876 const unsigned int global_num_cells,
3877 const std::string &file_basename) const;
3878
3887 void
3888 load_attached_data(const unsigned int global_first_cell,
3889 const unsigned int global_num_cells,
3890 const unsigned int local_num_cells,
3891 const std::string &file_basename,
3892 const unsigned int n_attached_deserialize_fixed,
3893 const unsigned int n_attached_deserialize_variable);
3894
3904 void
3906
3912 std::vector<typename internal::CellAttachedDataSerializer<dim, spacedim>::
3913 cell_relation_t>
3915
3921public:
3932 DeclException2(ExcInvalidLevel,
3933 int,
3934 int,
3935 << "You are requesting information from refinement level "
3936 << arg1
3937 << " of a triangulation, but this triangulation only has "
3938 << arg2 << " refinement levels. The given level " << arg1
3939 << " must be *less* than " << arg2 << '.');
3947 ExcTriangulationNotEmpty,
3948 int,
3949 int,
3950 << "You are trying to perform an operation on a triangulation "
3951 << "that is only allowed if the triangulation is currently empty. "
3952 << "However, it currently stores " << arg1 << " vertices and has "
3953 << "cells on " << arg2 << " levels.");
3959 DeclException0(ExcGridReadError);
3964 DeclException0(ExcFacesHaveNoLevel);
3970 DeclException1(ExcEmptyLevel,
3971 int,
3972 << "You tried to do something on level " << arg1
3973 << ", but this level is empty.");
3974
3982 DeclException1(ExcBoundaryIdNotFound,
3984 << "The given boundary_id " << arg1
3985 << " is not defined in this Triangulation!");
3986
3993 ExcInconsistentCoarseningFlags,
3994 "A cell is flagged for coarsening, but either not all of its siblings "
3995 "are active or flagged for coarsening as well. Please clean up all "
3996 "coarsen flags on your triangulation via "
3997 "Triangulation::prepare_coarsening_and_refinement() beforehand!");
3998
4001protected:
4007
4012 std::vector<ReferenceCell> reference_cells;
4013
4027 static void
4028 write_bool_vector(const unsigned int magic_number1,
4029 const std::vector<bool> &v,
4030 const unsigned int magic_number2,
4031 std::ostream &out);
4032
4037 static void
4038 read_bool_vector(const unsigned int magic_number1,
4039 std::vector<bool> &v,
4040 const unsigned int magic_number2,
4041 std::istream &in);
4042
4047 void
4049
4053 virtual void
4055
4056
4057private:
4062 std::unique_ptr<
4065
4072 std::vector<GridTools::PeriodicFacePair<cell_iterator>>
4074
4079 std::map<std::pair<cell_iterator, unsigned int>,
4080 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>>
4082
4093 TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
4096 using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
4097 using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
4098 using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
4099
4110 begin_raw(const unsigned int level = 0) const;
4111
4117 end_raw(const unsigned int level) const;
4118
4133 begin_raw_line(const unsigned int level = 0) const;
4134
4153 begin_line(const unsigned int level = 0) const;
4154
4173 begin_active_line(const unsigned int level = 0) const;
4174
4180 end_line() const;
4181
4209 begin_raw_quad(const unsigned int level = 0) const;
4210
4229 begin_quad(const unsigned int level = 0) const;
4230
4249 begin_active_quad(const unsigned int level = 0) const;
4250
4256 end_quad() const;
4257
4284 begin_raw_hex(const unsigned int level = 0) const;
4285
4304 begin_hex(const unsigned int level = 0) const;
4305
4324 begin_active_hex(const unsigned int level = 0) const;
4325
4331 end_hex() const;
4332
4349 void
4351
4355 void
4357
4364 void
4366
4370 void
4372
4376 void
4378
4394
4401 void
4403
4408 void
4410
4426 virtual unsigned int
4428 const types::coarse_cell_id coarse_cell_id) const;
4429
4430
4443 virtual types::coarse_cell_id
4445 const unsigned int coarse_cell_index) const;
4446
4451 std::vector<
4452 std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
4454
4460 std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
4462
4463
4467 std::vector<Point<spacedim>> vertices;
4468
4472 std::vector<bool> vertices_used;
4473
4478 std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4480
4485
4486
4492
4503
4518 std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4520
4521
4541 std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4543
4544 // make a couple of classes friends
4545 template <int, int, int>
4546 friend class TriaAccessorBase;
4547 template <int, int, int>
4548 friend class TriaAccessor;
4549 friend class TriaAccessor<0, 1, spacedim>;
4550
4551 friend class CellAccessor<dim, spacedim>;
4552
4553 friend struct ::internal::TriaAccessorImplementation::Implementation;
4554
4555 friend struct ::internal::TriangulationImplementation::Implementation;
4556 friend struct ::internal::TriangulationImplementation::
4557 ImplementationMixedMesh;
4558
4559 friend class ::internal::TriangulationImplementation::TriaObjects;
4560
4561 // explicitly check for sensible template arguments, but not on windows
4562 // because MSVC creates bogus warnings during normal compilation
4563#ifndef DEAL_II_MSVC
4564 static_assert(dim <= spacedim,
4565 "The dimension <dim> of a Triangulation must be less than or "
4566 "equal to the space dimension <spacedim> in which it lives.");
4567#endif
4568};
4569
4570
4571#ifndef DOXYGEN
4572
4573
4574
4575namespace internal
4576{
4577 namespace TriangulationImplementation
4578 {
4579 template <class Archive>
4580 void
4581 NumberCache<1>::serialize(Archive &ar, const unsigned int)
4582 {
4583 ar &n_levels;
4584 ar &n_lines &n_lines_level;
4585 ar &n_active_lines &n_active_lines_level;
4586 }
4587
4588
4589 template <class Archive>
4590 void
4591 NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4592 {
4593 this->NumberCache<1>::serialize(ar, version);
4594
4595 ar &n_quads &n_quads_level;
4596 ar &n_active_quads &n_active_quads_level;
4597 }
4598
4599
4600 template <class Archive>
4601 void
4602 NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4603 {
4604 this->NumberCache<2>::serialize(ar, version);
4605
4606 ar &n_hexes &n_hexes_level;
4607 ar &n_active_hexes &n_active_hexes_level;
4608 }
4609
4610 } // namespace TriangulationImplementation
4611} // namespace internal
4612
4613
4614template <int dim, int spacedim>
4617 const unsigned int index) const
4618{
4619 AssertIndexRange(index, vertices_used.size());
4620 return vertices_used[index];
4621}
4622
4623
4624
4625template <int dim, int spacedim>
4627inline unsigned int Triangulation<dim, spacedim>::n_levels() const
4628{
4629 return number_cache.n_levels;
4630}
4631
4632template <int dim, int spacedim>
4634inline unsigned int Triangulation<dim, spacedim>::n_global_levels() const
4635{
4636 return number_cache.n_levels;
4637}
4638
4639
4640template <int dim, int spacedim>
4642inline unsigned int Triangulation<dim, spacedim>::n_vertices() const
4643{
4644 return vertices.size();
4645}
4646
4647
4648
4649template <int dim, int spacedim>
4651inline const std::vector<Point<spacedim>>
4653{
4654 return vertices;
4655}
4656
4657
4658template <int dim, int spacedim>
4660template <class Archive>
4661void Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4662{
4663 // as discussed in the documentation, do not store the signals as
4664 // well as boundary and manifold description but everything else
4665 ar &smooth_grid;
4666
4667 unsigned int n_levels = levels.size();
4668 ar &n_levels;
4669 for (const auto &level : levels)
4670 ar &level;
4671
4672 // boost dereferences a nullptr when serializing a nullptr
4673 // at least up to 1.65.1. This causes problems with clang-5.
4674 // Therefore, work around it.
4675 bool faces_is_nullptr = (faces.get() == nullptr);
4676 ar &faces_is_nullptr;
4677 if (!faces_is_nullptr)
4678 ar &faces;
4679
4680 ar &vertices;
4681 ar &vertices_used;
4682
4683 ar &anisotropic_refinement;
4684 ar &number_cache;
4685
4686 ar &check_for_distorted_cells;
4687
4688 if (dim == 1)
4689 {
4690 ar &vertex_to_boundary_id_map_1d;
4691 ar &vertex_to_manifold_id_map_1d;
4692 }
4693}
4694
4695
4696
4697template <int dim, int spacedim>
4699template <class Archive>
4700void Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4701{
4702 // clear previous content. this also calls the respective signal
4703 clear();
4704
4705 // as discussed in the documentation, do not store the signals as
4706 // well as boundary and manifold description but everything else
4707 ar &smooth_grid;
4708
4709 unsigned int size;
4710 ar &size;
4711 levels.resize(size);
4712 for (auto &level_ : levels)
4713 {
4714 std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4715 ar &level;
4716 level_ = std::move(level);
4717 }
4718
4719 // Workaround for nullptr, see in save().
4720 bool faces_is_nullptr = true;
4721 ar &faces_is_nullptr;
4722 if (!faces_is_nullptr)
4723 ar &faces;
4724
4725 ar &vertices;
4726 ar &vertices_used;
4727
4728 ar &anisotropic_refinement;
4729 ar &number_cache;
4730
4731 // the levels do not serialize the active_cell_indices because
4732 // they are easy enough to rebuild upon re-loading data. do
4733 // this here. don't forget to first resize the fields appropriately
4734 {
4735 for (const auto &level : levels)
4736 {
4737 level->active_cell_indices.resize(level->refine_flags.size());
4738 level->global_active_cell_indices.resize(level->refine_flags.size());
4739 level->global_level_cell_indices.resize(level->refine_flags.size());
4740 }
4741 reset_cell_vertex_indices_cache();
4742 reset_active_cell_indices();
4743 reset_global_cell_indices();
4744 }
4745
4746 reset_policy();
4747
4748 bool my_check_for_distorted_cells;
4749 ar &my_check_for_distorted_cells;
4750
4751 Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4752 ExcMessage("The triangulation loaded into here must have the "
4753 "same setting with regard to reporting distorted "
4754 "cell as the one previously stored."));
4755
4756 if (dim == 1)
4757 {
4758 ar &vertex_to_boundary_id_map_1d;
4759 ar &vertex_to_manifold_id_map_1d;
4760 }
4761
4762 // trigger the create signal to indicate
4763 // that new content has been imported into
4764 // the triangulation
4765 signals.create();
4766}
4767
4768
4769
4770template <int dim, int spacedim>
4772inline unsigned int Triangulation<dim, spacedim>::
4774 const types::coarse_cell_id coarse_cell_id) const
4775{
4776 return coarse_cell_id;
4777}
4778
4779
4780
4781template <int dim, int spacedim>
4785 const unsigned int coarse_cell_index) const
4786{
4787 return coarse_cell_index;
4788}
4789
4790
4791
4792/* -------------- declaration of explicit specializations ------------- */
4793
4794template <>
4795unsigned int
4797template <>
4798unsigned int
4799Triangulation<1, 1>::n_quads(const unsigned int level) const;
4800template <>
4801unsigned int
4802Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4803template <>
4804unsigned int
4805Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4806template <>
4807unsigned int
4808Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4809template <>
4810unsigned int
4812template <>
4813unsigned int
4814Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4815template <>
4816unsigned int
4818template <>
4819unsigned int
4820Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4821template <>
4822unsigned int
4823Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4824template <>
4825unsigned int
4827template <>
4828unsigned int
4830template <>
4831unsigned int
4832Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4833template <>
4834unsigned int
4835Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4836
4837template <>
4838unsigned int
4840
4841
4842// -------------------------------------------------------------------
4843// -- Explicit specializations for codimension one grids
4844
4845
4846template <>
4847unsigned int
4849template <>
4850unsigned int
4851Triangulation<1, 2>::n_quads(const unsigned int level) const;
4852template <>
4853unsigned int
4854Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4855template <>
4856unsigned int
4857Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4858template <>
4859unsigned int
4860Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4861template <>
4862unsigned int
4863Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4864template <>
4865unsigned int
4867template <>
4868unsigned int
4870
4871// -------------------------------------------------------------------
4872// -- Explicit specializations for codimension two grids
4873
4874template <>
4875unsigned int
4877template <>
4878unsigned int
4879Triangulation<1, 3>::n_quads(const unsigned int level) const;
4880template <>
4881unsigned int
4882Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4883template <>
4884unsigned int
4885Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4886template <>
4887unsigned int
4888Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4889template <>
4890unsigned int
4891Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4892template <>
4893unsigned int
4895template <>
4896unsigned int
4898
4899template <>
4900bool
4902template <>
4903bool
4905template <>
4906bool
4908
4909
4910extern template class Triangulation<1, 1>;
4911extern template class Triangulation<1, 2>;
4912extern template class Triangulation<1, 3>;
4913extern template class Triangulation<2, 2>;
4914extern template class Triangulation<2, 3>;
4915extern template class Triangulation<3, 3>;
4916
4917#endif // DOXYGEN
4918
4920
4921// Include tria_accessor.h here, so that it is possible for an end
4922// user to use the iterators of Triangulation<dim> directly without
4923// the need to include tria_accessor.h separately. (Otherwise the
4924// iterators are an 'opaque' or 'incomplete' type.)
4926
4927#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:4096
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:4472
bool anisotropic_refinement
Definition tria.h:4484
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, unsigned char > > periodic_face_map
Definition tria.h:4081
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:4542
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:3914
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:4467
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:4491
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:4519
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:4012
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:4461
::CellStatus CellStatus
Definition tria.h:2228
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:4097
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:4502
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)
cell_iterator end() const
virtual bool has_hanging_nodes() const
std::vector< GridTools::PeriodicFacePair< cell_iterator > > periodic_face_pairs_level_0
Definition tria.h:4073
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:3864
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
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:4098
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
Definition tria.h:4479
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition tria.h:4006
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition tria.h:4064
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition tria.h:2513
internal::CellAttachedDataSerializer< dim, spacedim > data_serializer
Definition tria.h:3916
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:4453
unsigned int n_raw_hexs(const unsigned int level) const
unsigned int n_lines(const unsigned int level) const
unsigned int n_active_hexs() const
virtual unsigned int coarse_cell_id_to_coarse_cell_index(const types::coarse_cell_id coarse_cell_id) const
unsigned int n_quads(const unsigned int level) const
void load_user_flags(std::istream &in)
void reset_policy()
void save_coarsen_flags(std::ostream &out) const
active_face_iterator begin_active_face() const
void clear_user_flags_line()
raw_line_iterator begin_raw_line(const unsigned int level=0) const
static void write_bool_vector(const unsigned int magic_number1, const std::vector< bool > &v, const unsigned int magic_number2, std::ostream &out)
active_cell_iterator begin_active(const unsigned int level=0) const
void execute_coarsening()
void save_refine_flags(std::vector< bool > &v) const
std::vector< char > dest_data_variable
Definition tria.h:535
std::vector< char > src_data_fixed
Definition tria.h:525
std::vector< char > dest_data_fixed
Definition tria.h:526
std::vector< int > src_sizes_variable
Definition tria.h:532
std::vector< unsigned int > sizes_fixed_cumulative
Definition tria.h:519
typename std::pair< cell_iterator, CellStatus > cell_relation_t
Definition tria.h:393
std::vector< int > dest_sizes_variable
Definition tria.h:533
std::vector< char > src_data_variable
Definition tria.h:534
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
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:1691
typename IteratorSelector::active_quad_iterator active_quad_iterator
Definition tria.h:1682
typename IteratorSelector::active_hex_iterator active_hex_iterator
Definition tria.h:1702
typename IteratorSelector::quad_iterator quad_iterator
Definition tria.h:1667
typename IteratorSelector::line_iterator line_iterator
Definition tria.h:1643
typename IteratorSelector::active_line_iterator active_line_iterator
Definition tria.h:1658
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:2271
virtual ~DistortedCellList() noexcept override
boost::signals2::signal< void()> mesh_movement
Definition tria.h:2328
boost::signals2::signal< void()> post_distributed_load
Definition tria.h:2502
boost::signals2::signal< void()> pre_distributed_save
Definition tria.h:2480
boost::signals2::signal< unsigned int(const cell_iterator &, const ::CellStatus), CellWeightSum< unsigned int > > weight
Definition tria.h:2426
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition tria.h:2359
boost::signals2::signal< void()> pre_distributed_refinement
Definition tria.h:2439
boost::signals2::signal< void()> post_distributed_refinement
Definition tria.h:2458
boost::signals2::signal< void()> pre_refinement
Definition tria.h:2304
boost::signals2::signal< void()> any_change
Definition tria.h:2385
boost::signals2::signal< void()> create
Definition tria.h:2295
boost::signals2::signal< void()> clear
Definition tria.h:2374
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition tria.h:2349
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition tria.h:2339
boost::signals2::signal< void()> post_refinement
Definition tria.h:2311
boost::signals2::signal< void()> pre_partition
Definition tria.h:2319
boost::signals2::signal< void()> pre_distributed_repartition
Definition tria.h:2465
boost::signals2::signal< void()> pre_distributed_load
Definition tria.h:2495
boost::signals2::signal< void()> post_p4est_refinement
Definition tria.h:2448
boost::signals2::signal< void()> post_distributed_repartition
Definition tria.h:2472
boost::signals2::signal< void()> post_distributed_save
Definition tria.h:2487
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:354
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