Reference documentation for deal.II version 9.5.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// Copyright (C) 1998 - 2023 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_tria_h
17#define dealii_tria_h
18
19
20#include <deal.II/base/config.h>
21
25#include <deal.II/base/point.h>
28
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} // namespace internal
332
333
334/*------------------------------------------------------------------------*/
335
336
1134template <int dim, int spacedim = dim>
1137{
1138private:
1145
1146public:
1152 {
1157 none = 0x0,
1201 limit_level_difference_at_vertices = 0x1,
1222 eliminate_unrefined_islands = 0x2,
1238 patch_level_1 = 0x4,
1262 coarsest_level_1 = 0x8,
1287 allow_anisotropic_smoothing = 0x10,
1320 eliminate_refined_inner_islands = 0x100,
1325 eliminate_refined_boundary_islands = 0x200,
1331 do_not_produce_unrefined_islands = 0x400,
1332
1337 smoothing_on_refinement =
1338 (limit_level_difference_at_vertices | eliminate_unrefined_islands),
1343 smoothing_on_coarsening =
1344 (eliminate_refined_inner_islands | eliminate_refined_boundary_islands |
1345 do_not_produce_unrefined_islands),
1346
1352 maximum_smoothing = 0xffff ^ allow_anisotropic_smoothing
1354
1371
1377
1395
1409 using face_iterator = TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1410
1423 TriaActiveIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1424
1434
1449
1457 using line_iterator = typename IteratorSelector::line_iterator;
1458
1472 using active_line_iterator = typename IteratorSelector::active_line_iterator;
1473
1481 using quad_iterator = typename IteratorSelector::quad_iterator;
1482
1496 using active_quad_iterator = typename IteratorSelector::active_quad_iterator;
1497
1505 using hex_iterator = typename IteratorSelector::hex_iterator;
1506
1516 using active_hex_iterator = typename IteratorSelector::active_hex_iterator;
1517
1537 {
1544 virtual ~DistortedCellList() noexcept override;
1545
1550 std::list<typename Triangulation<dim, spacedim>::cell_iterator>
1551 distorted_cells;
1552 };
1553
1557 static constexpr unsigned int dimension = dim;
1558
1562 static constexpr unsigned int space_dimension = spacedim;
1563
1578 Triangulation(const MeshSmoothing smooth_grid = none,
1579 const bool check_for_distorted_cells = false);
1580
1596 Triangulation(const Triangulation<dim, spacedim> &) = delete;
1597
1604 Triangulation(Triangulation<dim, spacedim> &&tria) noexcept;
1605
1610 operator=(Triangulation<dim, spacedim> &&tria) noexcept;
1611
1615 virtual ~Triangulation() override;
1616
1623 virtual void
1624 clear();
1625
1630 virtual MPI_Comm
1631 get_communicator() const;
1632
1638 virtual const std::weak_ptr<const Utilities::MPI::Partitioner>
1639 global_active_cell_index_partitioner() const;
1640
1646 virtual const std::weak_ptr<const Utilities::MPI::Partitioner>
1647 global_level_cell_index_partitioner(const unsigned int level) const;
1648
1654 virtual void
1655 set_mesh_smoothing(const MeshSmoothing mesh_smoothing);
1656
1660 virtual const MeshSmoothing &
1661 get_mesh_smoothing() const;
1662
1685 void
1686 set_manifold(const types::manifold_id number,
1687 const Manifold<dim, spacedim> &manifold_object);
1688
1701 void
1702 reset_manifold(const types::manifold_id manifold_number);
1703
1715 void
1716 reset_all_manifolds();
1717
1726 void
1727 set_all_manifold_ids(const types::manifold_id number);
1728
1737 void
1738 set_all_manifold_ids_on_boundary(const types::manifold_id number);
1739
1749 void
1750 set_all_manifold_ids_on_boundary(const types::boundary_id b_id,
1751 const types::manifold_id number);
1752
1764 const Manifold<dim, spacedim> &
1765 get_manifold(const types::manifold_id number) const;
1766
1779 virtual std::vector<types::boundary_id>
1780 get_boundary_ids() const;
1781
1794 virtual std::vector<types::manifold_id>
1795 get_manifold_ids() const;
1796
1823 virtual void
1824 copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
1825
1873 virtual void
1874 create_triangulation(const std::vector<Point<spacedim>> &vertices,
1875 const std::vector<CellData<dim>> & cells,
1876 const SubCellData & subcelldata);
1877
1890 virtual void
1891 create_triangulation(
1892 const TriangulationDescription::Description<dim, spacedim>
1893 &construction_data);
1894
1904 virtual void
1905 create_triangulation_compatibility(
1906 const std::vector<Point<spacedim>> &vertices,
1907 const std::vector<CellData<dim>> & cells,
1908 const SubCellData & subcelldata);
1909
1916 void
1917 flip_all_direction_flags();
1918
1930 void
1931 set_all_refine_flags();
1932
1955 void
1956 refine_global(const unsigned int times = 1);
1957
1970 void
1971 coarsen_global(const unsigned int times = 1);
1972
2004 virtual void
2005 execute_coarsening_and_refinement();
2006
2036 virtual bool
2037 prepare_coarsening_and_refinement();
2038
2039 /*
2040 * @}
2041 */
2042
2057 {
2074 CELL_INVALID
2076
2082 template <typename T>
2084 {
2085 using result_type = T;
2086
2087 template <typename InputIterator>
2088 T
2089 operator()(InputIterator first, InputIterator last) const
2090 {
2091 return std::accumulate(first, last, T());
2092 }
2093 };
2094
2104 struct Signals
2105 {
2113 boost::signals2::signal<void()> create;
2114
2122 boost::signals2::signal<void()> pre_refinement;
2123
2129 boost::signals2::signal<void()> post_refinement;
2130
2137 boost::signals2::signal<void()> pre_partition;
2138
2146 boost::signals2::signal<void()> mesh_movement;
2147
2155 boost::signals2::signal<void(
2156 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2158
2165 boost::signals2::signal<void(
2166 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2168
2175 boost::signals2::signal<void(
2176 const Triangulation<dim, spacedim> &destination_tria)>
2178
2192 boost::signals2::signal<void()> clear;
2193
2203 boost::signals2::signal<void()> any_change;
2204
2240 boost::signals2::signal<unsigned int(const cell_iterator &,
2241 const CellStatus),
2244
2251 : cell_weight(weight)
2252 {}
2253
2258 {
2259 public:
2260 using signature_type = unsigned int(const cell_iterator &,
2261 const CellStatus);
2263
2264 using slot_function_type = boost::function<signature_type>;
2266 boost::signals2::slot<signature_type, slot_function_type>;
2267
2272 boost::signals2::signal<signature_type, combiner_type> &new_signal)
2273 : new_signal(new_signal)
2274 {}
2275
2280 {
2281 base_weight.disconnect();
2282 }
2283
2291 boost::signals2::connection
2293 const slot_type & slot,
2294 boost::signals2::connect_position position = boost::signals2::at_back)
2295 {
2296 if (base_weight.connected() == false)
2297 {
2298 base_weight = new_signal.connect(
2299 [](const cell_iterator &, const CellStatus) -> unsigned int {
2300 return 1000;
2301 });
2302 Assert(base_weight.connected() && new_signal.num_slots() == 1,
2304 }
2305
2306 return new_signal.connect(slot, position);
2307 }
2308
2314 std::size_t
2316 {
2317 return new_signal.num_slots() -
2318 static_cast<std::size_t>(base_weight.connected());
2319 }
2320
2325 bool
2326 empty() const
2327 {
2328 if (num_slots() == 0)
2329 {
2330 Assert(new_signal.num_slots() == 0, ExcInternalError());
2331 return true;
2332 }
2333 return false;
2334 }
2335
2342 template <typename S>
2344 disconnect(const S &connection)
2345 {
2346 new_signal.disconnect(connection);
2347
2348 if (num_slots() == 0)
2349 {
2350 Assert(base_weight.connected() && new_signal.num_slots() == 1,
2352 new_signal.disconnect(base_weight);
2353 }
2354 }
2355
2360 unsigned int
2361 operator()(const cell_iterator &iterator, const CellStatus status)
2362 {
2363 return new_signal(iterator, status);
2364 }
2365
2366 private:
2370 boost::signals2::connection base_weight;
2371
2375 boost::signals2::signal<signature_type, combiner_type> &new_signal;
2376 };
2377
2398
2410 boost::signals2::signal<void()> pre_distributed_refinement;
2411
2419 boost::signals2::signal<void()> post_p4est_refinement;
2420
2429 boost::signals2::signal<void()> post_distributed_refinement;
2430
2436 boost::signals2::signal<void()> pre_distributed_repartition;
2437
2443 boost::signals2::signal<void()> post_distributed_repartition;
2444
2451 boost::signals2::signal<void()> pre_distributed_save;
2452
2458 boost::signals2::signal<void()> post_distributed_save;
2459
2466 boost::signals2::signal<void()> pre_distributed_load;
2467
2473 boost::signals2::signal<void()> post_distributed_load;
2474 };
2475
2480
2481 /*
2482 * @}
2483 */
2484
2494 void
2495 save_refine_flags(std::ostream &out) const;
2496
2500 void
2501 save_refine_flags(std::vector<bool> &v) const;
2502
2506 void
2507 load_refine_flags(std::istream &in);
2508
2512 void
2513 load_refine_flags(const std::vector<bool> &v);
2514
2518 void
2519 save_coarsen_flags(std::ostream &out) const;
2520
2524 void
2525 save_coarsen_flags(std::vector<bool> &v) const;
2526
2530 void
2531 load_coarsen_flags(std::istream &out);
2532
2536 void
2537 load_coarsen_flags(const std::vector<bool> &v);
2538
2543 bool
2545
2546 /*
2547 * @}
2548 */
2549
2559 void
2561
2567 void
2568 save_user_flags(std::ostream &out) const;
2569
2575 void
2576 save_user_flags(std::vector<bool> &v) const;
2577
2582 void
2583 load_user_flags(std::istream &in);
2584
2589 void
2590 load_user_flags(const std::vector<bool> &v);
2591
2596 void
2598
2603 void
2604 save_user_flags_line(std::ostream &out) const;
2605
2611 void
2612 save_user_flags_line(std::vector<bool> &v) const;
2613
2618 void
2619 load_user_flags_line(std::istream &in);
2620
2625 void
2626 load_user_flags_line(const std::vector<bool> &v);
2627
2632 void
2634
2639 void
2640 save_user_flags_quad(std::ostream &out) const;
2641
2647 void
2648 save_user_flags_quad(std::vector<bool> &v) const;
2649
2654 void
2655 load_user_flags_quad(std::istream &in);
2656
2661 void
2662 load_user_flags_quad(const std::vector<bool> &v);
2663
2664
2669 void
2671
2676 void
2677 save_user_flags_hex(std::ostream &out) const;
2678
2684 void
2685 save_user_flags_hex(std::vector<bool> &v) const;
2686
2691 void
2692 load_user_flags_hex(std::istream &in);
2693
2698 void
2699 load_user_flags_hex(const std::vector<bool> &v);
2700
2706 void
2708
2714 void
2715 save_user_indices(std::vector<unsigned int> &v) const;
2716
2721 void
2722 load_user_indices(const std::vector<unsigned int> &v);
2723
2729 void
2730 save_user_pointers(std::vector<void *> &v) const;
2731
2736 void
2737 load_user_pointers(const std::vector<void *> &v);
2738
2744 void
2745 save_user_indices_line(std::vector<unsigned int> &v) const;
2746
2751 void
2752 load_user_indices_line(const std::vector<unsigned int> &v);
2753
2759 void
2760 save_user_indices_quad(std::vector<unsigned int> &v) const;
2761
2766 void
2767 load_user_indices_quad(const std::vector<unsigned int> &v);
2768
2774 void
2775 save_user_indices_hex(std::vector<unsigned int> &v) const;
2776
2781 void
2782 load_user_indices_hex(const std::vector<unsigned int> &v);
2788 void
2789 save_user_pointers_line(std::vector<void *> &v) const;
2790
2795 void
2796 load_user_pointers_line(const std::vector<void *> &v);
2797
2803 void
2804 save_user_pointers_quad(std::vector<void *> &v) const;
2805
2810 void
2811 load_user_pointers_quad(const std::vector<void *> &v);
2812
2818 void
2819 save_user_pointers_hex(std::vector<void *> &v) const;
2820
2825 void
2826 load_user_pointers_hex(const std::vector<void *> &v);
2827
2828 /*
2829 * @}
2830 */
2831
2855 begin(const unsigned int level = 0) const;
2856
2888 begin_active(const unsigned int level = 0) const;
2889
2895 end() const;
2896
2916 end(const unsigned int level) const;
2917
2938 end_active(const unsigned int level) const;
2939
2940
2945 last() const;
2946
2952
2966 create_cell_iterator(const CellId &cell_id) const;
2967
2984
3022
3039 cell_iterators_on_level(const unsigned int level) const;
3040
3057 active_cell_iterators_on_level(const unsigned int level) const;
3058
3059 /*
3060 * @}
3061 */
3062
3063 /*-------------------------------------------------------------------------*/
3064
3074 begin_face() const;
3075
3081
3087 end_face() const;
3088
3109
3110 /*
3111 * @}
3112 */
3113
3114 /*-------------------------------------------------------------------------*/
3115
3127
3135
3142 end_vertex() const;
3143
3144 /*
3145 * @}
3146 */
3147
3165 unsigned int
3166 n_lines() const;
3167
3171 unsigned int
3172 n_lines(const unsigned int level) const;
3173
3177 unsigned int
3179
3183 unsigned int
3184 n_active_lines(const unsigned int level) const;
3185
3189 unsigned int
3190 n_quads() const;
3191
3195 unsigned int
3196 n_quads(const unsigned int level) const;
3197
3201 unsigned int
3203
3207 unsigned int
3208 n_active_quads(const unsigned int level) const;
3209
3213 unsigned int
3214 n_hexs() const;
3215
3220 unsigned int
3221 n_hexs(const unsigned int level) const;
3222
3226 unsigned int
3228
3233 unsigned int
3234 n_active_hexs(const unsigned int level) const;
3235
3240 unsigned int
3241 n_cells() const;
3242
3247 unsigned int
3248 n_cells(const unsigned int level) const;
3249
3254 unsigned int
3256
3266
3267
3272 unsigned int
3273 n_active_cells(const unsigned int level) const;
3274
3279 virtual types::coarse_cell_id
3281
3287 unsigned int
3288 n_faces() const;
3289
3295 unsigned int
3297
3315 unsigned int
3316 n_levels() const;
3317
3324 virtual unsigned int
3326
3336 virtual bool
3338
3346 unsigned int
3347 n_vertices() const;
3348
3357 const std::vector<Point<spacedim>> &
3359
3364 unsigned int
3366
3370 bool
3371 vertex_used(const unsigned int index) const;
3372
3377 const std::vector<bool> &
3379
3391 unsigned int
3393
3400 virtual types::subdomain_id
3402
3414
3421
3422
3423 /*
3424 * @}
3425 */
3426
3441 unsigned int
3443
3453 unsigned int
3454 n_raw_lines(const unsigned int level) const;
3455
3465 unsigned int
3467
3477 unsigned int
3478 n_raw_quads(const unsigned int level) const;
3479
3489 unsigned int
3490 n_raw_hexs(const unsigned int level) const;
3491
3501 unsigned int
3502 n_raw_cells(const unsigned int level) const;
3503
3515 unsigned int
3517
3518 /*
3519 * @}
3520 */
3521
3530 virtual std::size_t
3532
3542 template <class Archive>
3543 void
3544 save(Archive &ar, const unsigned int version) const;
3545
3563 template <class Archive>
3564 void
3565 load(Archive &ar, const unsigned int version);
3566
3567
3583 virtual void
3585 const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3586
3590 const std::map<
3591 std::pair<cell_iterator, unsigned int>,
3592 std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3>>> &
3594
3599 const std::vector<ReferenceCell> &
3601
3606 bool
3608
3613 bool
3615
3621 bool
3623
3624#ifdef DOXYGEN
3630 template <class Archive>
3631 void
3632 serialize(Archive &archive, const unsigned int version);
3633#else
3634 // This macro defines the serialize() method that is compatible with
3635 // the templated save() and load() method that have been implemented.
3636 BOOST_SERIALIZATION_SPLIT_MEMBER()
3637#endif
3638
3649 DeclException2(ExcInvalidLevel,
3650 int,
3651 int,
3652 << "You are requesting information from refinement level "
3653 << arg1
3654 << " of a triangulation, but this triangulation only has "
3655 << arg2 << " refinement levels. The given level " << arg1
3656 << " must be *less* than " << arg2 << '.');
3664 ExcTriangulationNotEmpty,
3665 int,
3666 int,
3667 << "You are trying to perform an operation on a triangulation "
3668 << "that is only allowed if the triangulation is currently empty. "
3669 << "However, it currently stores " << arg1 << " vertices and has "
3670 << "cells on " << arg2 << " levels.");
3676 DeclException0(ExcGridReadError);
3681 DeclException0(ExcFacesHaveNoLevel);
3687 DeclException1(ExcEmptyLevel,
3688 int,
3689 << "You tried to do something on level " << arg1
3690 << ", but this level is empty.");
3696 DeclException0(ExcNonOrientableTriangulation);
3697
3705 DeclException1(ExcBoundaryIdNotFound,
3707 << "The given boundary_id " << arg1
3708 << " is not defined in this Triangulation!");
3709
3716 ExcInconsistentCoarseningFlags,
3717 "A cell is flagged for coarsening, but either not all of its siblings "
3718 "are active or flagged for coarsening as well. Please clean up all "
3719 "coarsen flags on your triangulation via "
3720 "Triangulation::prepare_coarsening_and_refinement() beforehand!");
3721
3722 /*
3723 * @}
3724 */
3725
3726protected:
3732
3737 std::vector<ReferenceCell> reference_cells;
3738
3752 static void
3753 write_bool_vector(const unsigned int magic_number1,
3754 const std::vector<bool> &v,
3755 const unsigned int magic_number2,
3756 std::ostream & out);
3757
3762 static void
3763 read_bool_vector(const unsigned int magic_number1,
3764 std::vector<bool> &v,
3765 const unsigned int magic_number2,
3766 std::istream & in);
3767
3772 void
3774
3778 virtual void
3780
3781
3782private:
3787 std::unique_ptr<
3790
3797 std::vector<GridTools::PeriodicFacePair<cell_iterator>>
3799
3804 std::map<std::pair<cell_iterator, unsigned int>,
3805 std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3>>>
3807
3823 TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
3826 using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
3827 using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
3828 using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
3829
3835 begin_raw(const unsigned int level = 0) const;
3836
3842 end_raw(const unsigned int level) const;
3843
3844 /*
3845 * @}
3846 */
3847
3860 begin_raw_line(const unsigned int level = 0) const;
3861
3880 begin_line(const unsigned int level = 0) const;
3881
3900 begin_active_line(const unsigned int level = 0) const;
3901
3907 end_line() const;
3908
3909 /*
3910 * @}
3911 */
3912
3938 begin_raw_quad(const unsigned int level = 0) const;
3939
3958 begin_quad(const unsigned int level = 0) const;
3959
3978 begin_active_quad(const unsigned int level = 0) const;
3979
3985 end_quad() const;
3986
3987 /*
3988 * @}
3989 */
3990
4015 begin_raw_hex(const unsigned int level = 0) const;
4016
4035 begin_hex(const unsigned int level = 0) const;
4036
4055 begin_active_hex(const unsigned int level = 0) const;
4056
4062 end_hex() const;
4063
4064 /*
4065 * @}
4066 */
4067
4068
4082 void
4084
4088 void
4090
4097 void
4099
4103 void
4105
4109 void
4111
4127
4134 void
4136
4141 void
4143
4157 virtual unsigned int
4159 const types::coarse_cell_id coarse_cell_id) const;
4160
4161
4174 virtual types::coarse_cell_id
4176 const unsigned int coarse_cell_index) const;
4177
4182 std::vector<
4183 std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
4185
4191 std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
4193
4194
4198 std::vector<Point<spacedim>> vertices;
4199
4203 std::vector<bool> vertices_used;
4204
4209 std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4211
4216
4217
4223
4234
4249 std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4251
4252
4272 std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4274
4275 // make a couple of classes friends
4276 template <int, int, int>
4277 friend class TriaAccessorBase;
4278 template <int, int, int>
4279 friend class TriaAccessor;
4280 friend class TriaAccessor<0, 1, spacedim>;
4281
4282 friend class CellAccessor<dim, spacedim>;
4283
4284 friend struct ::internal::TriaAccessorImplementation::Implementation;
4285
4286 friend struct ::internal::TriangulationImplementation::Implementation;
4287 friend struct ::internal::TriangulationImplementation::
4288 ImplementationMixedMesh;
4289
4290 friend class ::internal::TriangulationImplementation::TriaObjects;
4291
4292 // explicitly check for sensible template arguments, but not on windows
4293 // because MSVC creates bogus warnings during normal compilation
4294#ifndef DEAL_II_MSVC
4295 static_assert(dim <= spacedim,
4296 "The dimension <dim> of a Triangulation must be less than or "
4297 "equal to the space dimension <spacedim> in which it lives.");
4298#endif
4299};
4300
4301
4302#ifndef DOXYGEN
4303
4304
4305
4306namespace internal
4307{
4308 namespace TriangulationImplementation
4309 {
4310 template <class Archive>
4311 void
4312 NumberCache<1>::serialize(Archive &ar, const unsigned int)
4313 {
4314 ar &n_levels;
4315 ar &n_lines &n_lines_level;
4316 ar &n_active_lines &n_active_lines_level;
4317 }
4318
4319
4320 template <class Archive>
4321 void
4322 NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4323 {
4324 this->NumberCache<1>::serialize(ar, version);
4325
4326 ar &n_quads &n_quads_level;
4327 ar &n_active_quads &n_active_quads_level;
4328 }
4329
4330
4331 template <class Archive>
4332 void
4333 NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4334 {
4335 this->NumberCache<2>::serialize(ar, version);
4336
4337 ar &n_hexes &n_hexes_level;
4338 ar &n_active_hexes &n_active_hexes_level;
4339 }
4340
4341 } // namespace TriangulationImplementation
4342} // namespace internal
4343
4344
4345template <int dim, int spacedim>
4348 const unsigned int index) const
4349{
4350 AssertIndexRange(index, vertices_used.size());
4351 return vertices_used[index];
4352}
4353
4354
4355
4356template <int dim, int spacedim>
4358inline unsigned int Triangulation<dim, spacedim>::n_levels() const
4359{
4360 return number_cache.n_levels;
4361}
4362
4363template <int dim, int spacedim>
4365inline unsigned int Triangulation<dim, spacedim>::n_global_levels() const
4366{
4367 return number_cache.n_levels;
4368}
4369
4370
4371template <int dim, int spacedim>
4373inline unsigned int Triangulation<dim, spacedim>::n_vertices() const
4374{
4375 return vertices.size();
4376}
4377
4378
4379
4380template <int dim, int spacedim>
4382inline const std::vector<Point<spacedim>>
4384{
4385 return vertices;
4386}
4387
4388
4389template <int dim, int spacedim>
4391template <class Archive>
4392void Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4393{
4394 // as discussed in the documentation, do not store the signals as
4395 // well as boundary and manifold description but everything else
4396 ar &smooth_grid;
4397
4398 unsigned int n_levels = levels.size();
4399 ar & n_levels;
4400 for (const auto &level : levels)
4401 ar &level;
4402
4403 // boost dereferences a nullptr when serializing a nullptr
4404 // at least up to 1.65.1. This causes problems with clang-5.
4405 // Therefore, work around it.
4406 bool faces_is_nullptr = (faces.get() == nullptr);
4407 ar & faces_is_nullptr;
4408 if (!faces_is_nullptr)
4409 ar &faces;
4410
4411 ar &vertices;
4412 ar &vertices_used;
4413
4414 ar &anisotropic_refinement;
4415 ar &number_cache;
4416
4417 ar &check_for_distorted_cells;
4418
4419 if (dim == 1)
4420 {
4421 ar &vertex_to_boundary_id_map_1d;
4422 ar &vertex_to_manifold_id_map_1d;
4423 }
4424}
4425
4426
4427
4428template <int dim, int spacedim>
4430template <class Archive>
4431void Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4432{
4433 // clear previous content. this also calls the respective signal
4434 clear();
4435
4436 // as discussed in the documentation, do not store the signals as
4437 // well as boundary and manifold description but everything else
4438 ar &smooth_grid;
4439
4440 unsigned int size;
4441 ar & size;
4442 levels.resize(size);
4443 for (auto &level_ : levels)
4444 {
4445 std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4446 ar & level;
4447 level_ = std::move(level);
4448 }
4449
4450 // Workaround for nullptr, see in save().
4451 bool faces_is_nullptr = true;
4452 ar & faces_is_nullptr;
4453 if (!faces_is_nullptr)
4454 ar &faces;
4455
4456 ar &vertices;
4457 ar &vertices_used;
4458
4459 ar &anisotropic_refinement;
4460 ar &number_cache;
4461
4462 // the levels do not serialize the active_cell_indices because
4463 // they are easy enough to rebuild upon re-loading data. do
4464 // this here. don't forget to first resize the fields appropriately
4465 {
4466 for (auto &level : levels)
4467 {
4468 level->active_cell_indices.resize(level->refine_flags.size());
4469 level->global_active_cell_indices.resize(level->refine_flags.size());
4470 level->global_level_cell_indices.resize(level->refine_flags.size());
4471 }
4472 reset_cell_vertex_indices_cache();
4473 reset_active_cell_indices();
4474 reset_global_cell_indices();
4475 }
4476
4477 reset_policy();
4478
4479 bool my_check_for_distorted_cells;
4480 ar & my_check_for_distorted_cells;
4481
4482 Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4483 ExcMessage("The triangulation loaded into here must have the "
4484 "same setting with regard to reporting distorted "
4485 "cell as the one previously stored."));
4486
4487 if (dim == 1)
4488 {
4489 ar &vertex_to_boundary_id_map_1d;
4490 ar &vertex_to_manifold_id_map_1d;
4491 }
4492
4493 // trigger the create signal to indicate
4494 // that new content has been imported into
4495 // the triangulation
4496 signals.create();
4497}
4498
4499
4500
4501template <int dim, int spacedim>
4503inline unsigned int Triangulation<dim, spacedim>::
4505 const types::coarse_cell_id coarse_cell_id) const
4506{
4507 return coarse_cell_id;
4508}
4509
4510
4511
4512template <int dim, int spacedim>
4516 const unsigned int coarse_cell_index) const
4517{
4518 return coarse_cell_index;
4519}
4520
4521
4522
4523/* -------------- declaration of explicit specializations ------------- */
4524
4525template <>
4526unsigned int
4528template <>
4529unsigned int
4530Triangulation<1, 1>::n_quads(const unsigned int level) const;
4531template <>
4532unsigned int
4533Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4534template <>
4535unsigned int
4536Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4537template <>
4538unsigned int
4539Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4540template <>
4541unsigned int
4543template <>
4544unsigned int
4545Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4546template <>
4547unsigned int
4549template <>
4550unsigned int
4551Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4552template <>
4553unsigned int
4554Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4555template <>
4556unsigned int
4558template <>
4559unsigned int
4561template <>
4562unsigned int
4563Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4564template <>
4565unsigned int
4566Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4567
4568template <>
4569unsigned int
4571
4572
4573// -------------------------------------------------------------------
4574// -- Explicit specializations for codimension one grids
4575
4576
4577template <>
4578unsigned int
4580template <>
4581unsigned int
4582Triangulation<1, 2>::n_quads(const unsigned int level) const;
4583template <>
4584unsigned int
4585Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4586template <>
4587unsigned int
4588Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4589template <>
4590unsigned int
4591Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4592template <>
4593unsigned int
4594Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4595template <>
4596unsigned int
4598template <>
4599unsigned int
4601
4602// -------------------------------------------------------------------
4603// -- Explicit specializations for codimension two grids
4604
4605template <>
4606unsigned int
4608template <>
4609unsigned int
4610Triangulation<1, 3>::n_quads(const unsigned int level) const;
4611template <>
4612unsigned int
4613Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4614template <>
4615unsigned int
4616Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4617template <>
4618unsigned int
4619Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4620template <>
4621unsigned int
4622Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4623template <>
4624unsigned int
4626template <>
4627unsigned int
4629
4630template <>
4631bool
4633template <>
4634bool
4636template <>
4637bool
4639
4640
4641extern template class Triangulation<1, 1>;
4642extern template class Triangulation<1, 2>;
4643extern template class Triangulation<1, 3>;
4644extern template class Triangulation<2, 2>;
4645extern template class Triangulation<2, 3>;
4646extern template class Triangulation<3, 3>;
4647
4648#endif // DOXYGEN
4649
4651
4652// Include tria_accessor.h here, so that it is possible for an end
4653// user to use the iterators of Triangulation<dim> directly without
4654// the need to include tria_accessor.h separately. (Otherwise the
4655// iterators are an 'opaque' or 'incomplete' type.)
4657
4658#endif
Definition point.h:112
unsigned int operator()(const cell_iterator &iterator, const CellStatus status)
Definition tria.h:2361
boost::signals2::slot< signature_type, slot_function_type > slot_type
Definition tria.h:2266
std::size_t num_slots() const
Definition tria.h:2315
boost::function< signature_type > slot_function_type
Definition tria.h:2264
boost::signals2::connection connect(const slot_type &slot, boost::signals2::connect_position position=boost::signals2::at_back)
Definition tria.h:2292
boost::signals2::signal< signature_type, combiner_type > & new_signal
Definition tria.h:2375
LegacySignal(boost::signals2::signal< signature_type, combiner_type > &new_signal)
Definition tria.h:2271
unsigned int(const cell_iterator &, const CellStatus) signature_type
Definition tria.h:2261
boost::signals2::connection base_weight
Definition tria.h:2370
void disconnect(const S &connection)
Definition tria.h:2344
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:3826
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:4203
bool anisotropic_refinement
Definition tria.h:4215
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:4273
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
const Triangulation< dim, spacedim > & get_triangulation() const
void reset_global_cell_indices()
face_iterator end_face() const
void reset_active_cell_indices()
cell_iterator create_cell_iterator(const CellId &cell_id) const
cell_iterator begin(const unsigned int level=0) const
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > periodic_face_map
Definition tria.h:3806
void fix_coarsen_flags()
unsigned int n_active_cells(const unsigned int level) const
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > & get_periodic_face_map() 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:4198
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:4222
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:4250
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:3737
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:4192
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:3827
void load_user_pointers(const std::vector< void * > &v)
::internal::TriangulationImplementation::NumberCache< dim > number_cache
Definition tria.h:4233
void load_user_flags_quad(const std::vector< bool > &v)
void load_user_flags_line(const std::vector< bool > &v)
void save_user_indices_hex(std::vector< unsigned int > &v) const
DistortedCellList execute_refinement()
active_line_iterator begin_active_line(const unsigned int level=0) const
void save_user_indices_quad(std::vector< unsigned int > &v) const
void load_user_pointers_hex(const std::vector< void * > &v)
cell_iterator end() const
virtual bool has_hanging_nodes() const
std::vector< GridTools::PeriodicFacePair< cell_iterator > > periodic_face_pairs_level_0
Definition tria.h:3798
unsigned int n_raw_cells(const unsigned int level) const
cell_iterator end(const unsigned int level) const
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const
unsigned int n_raw_quads(const unsigned int level) const
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:3828
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
Definition tria.h:4210
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition tria.h:3731
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition tria.h:3789
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition tria.h:2479
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 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:4184
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
#define DEAL_II_DEPRECATED
Definition config.h:172
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:160
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
Point< 3 > vertices[4]
Point< 2 > first
Definition grid_out.cc:4615
unsigned int level
Definition grid_out.cc:4618
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:465
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:533
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:488
static ::ExceptionBase & ExcInternalError()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:510
static ::ExceptionBase & ExcMessage(std::string arg1)
typename IteratorSelector::hex_iterator hex_iterator
Definition tria.h:1505
typename IteratorSelector::active_quad_iterator active_quad_iterator
Definition tria.h:1496
typename IteratorSelector::active_hex_iterator active_hex_iterator
Definition tria.h:1516
typename IteratorSelector::quad_iterator quad_iterator
Definition tria.h:1481
typename IteratorSelector::line_iterator line_iterator
Definition tria.h:1457
typename IteratorSelector::active_line_iterator active_line_iterator
Definition tria.h:1472
STL namespace.
Definition types.h:33
global_cell_index coarse_cell_id
Definition types.h:127
T operator()(InputIterator first, InputIterator last) const
Definition tria.h:2089
virtual ~DistortedCellList() noexcept override
boost::signals2::signal< unsigned int(const cell_iterator &, const CellStatus), CellWeightSum< unsigned int > > weight
Definition tria.h:2243
boost::signals2::signal< void()> mesh_movement
Definition tria.h:2146
boost::signals2::signal< void()> post_distributed_load
Definition tria.h:2473
boost::signals2::signal< void()> pre_distributed_save
Definition tria.h:2451
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition tria.h:2177
boost::signals2::signal< void()> pre_distributed_refinement
Definition tria.h:2410
boost::signals2::signal< void()> post_distributed_refinement
Definition tria.h:2429
boost::signals2::signal< void()> pre_refinement
Definition tria.h:2122
boost::signals2::signal< void()> any_change
Definition tria.h:2203
boost::signals2::signal< void()> create
Definition tria.h:2113
boost::signals2::signal< void()> clear
Definition tria.h:2192
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition tria.h:2167
LegacySignal cell_weight
Definition tria.h:2397
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition tria.h:2157
boost::signals2::signal< void()> post_refinement
Definition tria.h:2129
boost::signals2::signal< void()> pre_partition
Definition tria.h:2137
boost::signals2::signal< void()> pre_distributed_repartition
Definition tria.h:2436
boost::signals2::signal< void()> pre_distributed_load
Definition tria.h:2466
boost::signals2::signal< void()> post_p4est_refinement
Definition tria.h:2419
boost::signals2::signal< void()> post_distributed_repartition
Definition tria.h:2443
boost::signals2::signal< void()> post_distributed_save
Definition tria.h:2458
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
const ::Triangulation< dim, spacedim > & tria