deal.II version GIT relicensing-2055-gd703fda1ad 2024-10-30 20:00:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
tria.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_tria_h
16#define dealii_tria_h
17
18
19#include <deal.II/base/config.h>
20
26#include <deal.II/base/point.h>
27
33
34#include <boost/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_mpi_communicator() const;
1818
1825 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
1826 "Access the MPI communicator with get_mpi_communicator() instead.")
1827 MPI_Comm
1828 get_communicator() const;
1829
1835 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1836 global_active_cell_index_partitioner() const;
1837
1843 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1844 global_level_cell_index_partitioner(const unsigned int level) const;
1845
1850 virtual void
1851 set_mesh_smoothing(const MeshSmoothing mesh_smoothing);
1852
1856 virtual const MeshSmoothing &
1857 get_mesh_smoothing() const;
1858
1881 void
1882 set_manifold(const types::manifold_id number,
1883 const Manifold<dim, spacedim> &manifold_object);
1884
1901 void
1902 reset_manifold(const types::manifold_id manifold_number);
1903
1920 void
1921 reset_all_manifolds();
1922
1931 void
1932 set_all_manifold_ids(const types::manifold_id number);
1933
1942 void
1943 set_all_manifold_ids_on_boundary(const types::manifold_id number);
1944
1954 void
1955 set_all_manifold_ids_on_boundary(const types::boundary_id b_id,
1956 const types::manifold_id number);
1957
1969 const Manifold<dim, spacedim> &
1970 get_manifold(const types::manifold_id number) const;
1971
1982 virtual std::vector<types::boundary_id>
1983 get_boundary_ids() const;
1984
1997 virtual std::vector<types::manifold_id>
1998 get_manifold_ids() const;
1999
2021 virtual void
2022 copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
2023
2072 virtual void
2073 create_triangulation(const std::vector<Point<spacedim>> &vertices,
2074 const std::vector<CellData<dim>> &cells,
2075 const SubCellData &subcelldata);
2076
2089 virtual void
2090 create_triangulation(
2091 const TriangulationDescription::Description<dim, spacedim>
2092 &construction_data);
2093
2102 void
2103 flip_all_direction_flags();
2104
2116 void
2117 set_all_refine_flags();
2118
2141 void
2142 refine_global(const unsigned int times = 1);
2143
2156 void
2157 coarsen_global(const unsigned int times = 1);
2158
2190 virtual void
2191 execute_coarsening_and_refinement();
2192
2222 virtual bool
2223 prepare_coarsening_and_refinement();
2224
2240
2245 static constexpr auto CELL_PERSIST DEAL_II_DEPRECATED =
2247
2252 static constexpr auto CELL_REFINE DEAL_II_DEPRECATED =
2254
2259 static constexpr auto CELL_COARSEN DEAL_II_DEPRECATED =
2261
2266 static constexpr auto CELL_INVALID DEAL_II_DEPRECATED =
2268
2269
2275 template <typename T>
2277 {
2278 using result_type = T;
2279
2280 template <typename InputIterator>
2281 T
2282 operator()(InputIterator first, InputIterator last) const
2283 {
2284 return std::accumulate(first, last, T());
2285 }
2286 };
2287
2297 struct Signals
2298 {
2306 boost::signals2::signal<void()> create;
2307
2315 boost::signals2::signal<void()> pre_refinement;
2316
2322 boost::signals2::signal<void()> post_refinement;
2323
2330 boost::signals2::signal<void()> pre_partition;
2331
2339 boost::signals2::signal<void()> mesh_movement;
2340
2348 boost::signals2::signal<void(
2349 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2351
2358 boost::signals2::signal<void(
2359 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2361
2368 boost::signals2::signal<void(
2369 const Triangulation<dim, spacedim> &destination_tria)>
2371
2385 boost::signals2::signal<void()> clear;
2386
2396 boost::signals2::signal<void()> any_change;
2397
2434 boost::signals2::signal<unsigned int(const cell_iterator &,
2435 const ::CellStatus),
2438
2450 boost::signals2::signal<void()> pre_distributed_refinement;
2451
2459 boost::signals2::signal<void()> post_p4est_refinement;
2460
2469 boost::signals2::signal<void()> post_distributed_refinement;
2470
2476 boost::signals2::signal<void()> pre_distributed_repartition;
2477
2483 boost::signals2::signal<void()> post_distributed_repartition;
2484
2491 boost::signals2::signal<void()> pre_distributed_save;
2492
2498 boost::signals2::signal<void()> post_distributed_save;
2499
2506 boost::signals2::signal<void()> pre_distributed_load;
2507
2513 boost::signals2::signal<void()> post_distributed_load;
2514 };
2515
2525
2537 void
2538 save_refine_flags(std::ostream &out) const;
2539
2543 void
2544 save_refine_flags(std::vector<bool> &v) const;
2545
2549 void
2550 load_refine_flags(std::istream &in);
2551
2555 void
2556 load_refine_flags(const std::vector<bool> &v);
2557
2561 void
2562 save_coarsen_flags(std::ostream &out) const;
2563
2567 void
2568 save_coarsen_flags(std::vector<bool> &v) const;
2569
2573 void
2574 load_coarsen_flags(std::istream &out);
2575
2579 void
2580 load_coarsen_flags(const std::vector<bool> &v);
2581
2586 bool
2588
2600 void
2602
2608 void
2609 save_user_flags(std::ostream &out) const;
2610
2616 void
2617 save_user_flags(std::vector<bool> &v) const;
2618
2623 void
2624 load_user_flags(std::istream &in);
2625
2630 void
2631 load_user_flags(const std::vector<bool> &v);
2632
2637 void
2639
2644 void
2645 save_user_flags_line(std::ostream &out) const;
2646
2652 void
2653 save_user_flags_line(std::vector<bool> &v) const;
2654
2659 void
2660 load_user_flags_line(std::istream &in);
2661
2666 void
2667 load_user_flags_line(const std::vector<bool> &v);
2668
2673 void
2675
2680 void
2681 save_user_flags_quad(std::ostream &out) const;
2682
2688 void
2689 save_user_flags_quad(std::vector<bool> &v) const;
2690
2695 void
2696 load_user_flags_quad(std::istream &in);
2697
2702 void
2703 load_user_flags_quad(const std::vector<bool> &v);
2704
2705
2710 void
2712
2717 void
2718 save_user_flags_hex(std::ostream &out) const;
2719
2725 void
2726 save_user_flags_hex(std::vector<bool> &v) const;
2727
2732 void
2733 load_user_flags_hex(std::istream &in);
2734
2739 void
2740 load_user_flags_hex(const std::vector<bool> &v);
2741
2747 void
2749
2755 void
2756 save_user_indices(std::vector<unsigned int> &v) const;
2757
2762 void
2763 load_user_indices(const std::vector<unsigned int> &v);
2764
2770 void
2771 save_user_pointers(std::vector<void *> &v) const;
2772
2777 void
2778 load_user_pointers(const std::vector<void *> &v);
2779
2785 void
2786 save_user_indices_line(std::vector<unsigned int> &v) const;
2787
2792 void
2793 load_user_indices_line(const std::vector<unsigned int> &v);
2794
2800 void
2801 save_user_indices_quad(std::vector<unsigned int> &v) const;
2802
2807 void
2808 load_user_indices_quad(const std::vector<unsigned int> &v);
2809
2815 void
2816 save_user_indices_hex(std::vector<unsigned int> &v) const;
2817
2822 void
2823 load_user_indices_hex(const std::vector<unsigned int> &v);
2829 void
2830 save_user_pointers_line(std::vector<void *> &v) const;
2831
2836 void
2837 load_user_pointers_line(const std::vector<void *> &v);
2838
2844 void
2845 save_user_pointers_quad(std::vector<void *> &v) const;
2846
2851 void
2852 load_user_pointers_quad(const std::vector<void *> &v);
2853
2859 void
2860 save_user_pointers_hex(std::vector<void *> &v) const;
2861
2866 void
2867 load_user_pointers_hex(const std::vector<void *> &v);
2868
2894 begin(const unsigned int level = 0) const;
2895
2927 begin_active(const unsigned int level = 0) const;
2928
2934 end() const;
2935
2955 end(const unsigned int level) const;
2956
2977 end_active(const unsigned int level) const;
2978
2979
2984 last() const;
2985
2991
3000 create_cell_iterator(const CellId &cell_id) const;
3001
3012 bool
3013 contains_cell(const CellId &cell_id) const;
3033
3071
3088 cell_iterators_on_level(const unsigned int level) const;
3089
3106 active_cell_iterators_on_level(const unsigned int level) const;
3107
3110 /*-------------------------------------------------------------------------*/
3111
3121 begin_face() const;
3122
3128
3134 end_face() const;
3135
3156
3159 /*-------------------------------------------------------------------------*/
3160
3172
3180
3187 end_vertex() const;
3188
3208 unsigned int
3209 n_lines() const;
3210
3214 unsigned int
3215 n_lines(const unsigned int level) const;
3216
3220 unsigned int
3222
3226 unsigned int
3227 n_active_lines(const unsigned int level) const;
3228
3232 unsigned int
3233 n_quads() const;
3234
3238 unsigned int
3239 n_quads(const unsigned int level) const;
3240
3244 unsigned int
3246
3250 unsigned int
3251 n_active_quads(const unsigned int level) const;
3252
3256 unsigned int
3257 n_hexs() const;
3258
3263 unsigned int
3264 n_hexs(const unsigned int level) const;
3265
3269 unsigned int
3271
3276 unsigned int
3277 n_active_hexs(const unsigned int level) const;
3278
3283 unsigned int
3284 n_cells() const;
3285
3290 unsigned int
3291 n_cells(const unsigned int level) const;
3292
3297 unsigned int
3299
3309
3310
3315 unsigned int
3316 n_active_cells(const unsigned int level) const;
3317
3322 virtual types::coarse_cell_id
3324
3330 unsigned int
3331 n_faces() const;
3332
3338 unsigned int
3340
3358 unsigned int
3359 n_levels() const;
3360
3367 virtual unsigned int
3369
3379 virtual bool
3381
3389 unsigned int
3390 n_vertices() const;
3391
3400 const std::vector<Point<spacedim>> &
3402
3407 unsigned int
3409
3413 bool
3414 vertex_used(const unsigned int index) const;
3415
3420 const std::vector<bool> &
3422
3434 unsigned int
3436
3443 virtual types::subdomain_id
3445
3457
3464
3465
3482 unsigned int
3484
3494 unsigned int
3495 n_raw_lines(const unsigned int level) const;
3496
3506 unsigned int
3508
3518 unsigned int
3519 n_raw_quads(const unsigned int level) const;
3520
3530 unsigned int
3531 n_raw_hexs(const unsigned int level) const;
3532
3542 unsigned int
3543 n_raw_cells(const unsigned int level) const;
3544
3556 unsigned int
3558
3569 virtual std::size_t
3571
3581 template <class Archive>
3582 void
3583 save(Archive &ar, const unsigned int version) const;
3584
3602 template <class Archive>
3603 void
3604 load(Archive &ar, const unsigned int version);
3605
3606
3615 virtual void
3616 save(const std::string &file_basename) const;
3617
3621 virtual void
3622 load(const std::string &file_basename);
3623
3624
3640 virtual void
3642 const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3643
3647 const std::map<
3648 std::pair<cell_iterator, unsigned int>,
3649 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>> &
3651
3656 const std::vector<ReferenceCell> &
3658
3663 bool
3665
3670 bool
3672
3678 bool
3680
3681#ifdef DOXYGEN
3689 template <class Archive>
3690 void
3691 serialize(Archive &archive, const unsigned int version);
3692#else
3693 // This macro defines the serialize() method that is compatible with
3694 // the templated save() and load() method that have been implemented.
3695 BOOST_SERIALIZATION_SPLIT_MEMBER()
3696#endif
3697
3702public:
3810 unsigned int
3812 const std::function<std::vector<char>(const cell_iterator &,
3813 const ::CellStatus)>
3814 &pack_callback,
3815 const bool returns_variable_size_data);
3816
3866 void
3868 const unsigned int handle,
3869 const std::function<
3870 void(const cell_iterator &,
3871 const ::CellStatus,
3872 const boost::iterator_range<std::vector<char>::const_iterator> &)>
3873 &unpack_callback);
3874
3876
3877protected:
3885 void
3886 save_attached_data(const unsigned int global_first_cell,
3887 const unsigned int global_num_cells,
3888 const std::string &file_basename) const;
3889
3898 void
3899 load_attached_data(const unsigned int global_first_cell,
3900 const unsigned int global_num_cells,
3901 const unsigned int local_num_cells,
3902 const std::string &file_basename,
3903 const unsigned int n_attached_deserialize_fixed,
3904 const unsigned int n_attached_deserialize_variable);
3905
3915 void
3917
3923 void
3925
3926
3931 void
3933
3939 std::vector<typename internal::CellAttachedDataSerializer<dim, spacedim>::
3940 cell_relation_t>
3942
3948public:
3959 DeclException2(ExcInvalidLevel,
3960 int,
3961 int,
3962 << "You are requesting information from refinement level "
3963 << arg1
3964 << " of a triangulation, but this triangulation only has "
3965 << arg2 << " refinement levels. The given level " << arg1
3966 << " must be *less* than " << arg2 << '.');
3974 ExcTriangulationNotEmpty,
3975 int,
3976 int,
3977 << "You are trying to perform an operation on a triangulation "
3978 << "that is only allowed if the triangulation is currently empty. "
3979 << "However, it currently stores " << arg1 << " vertices and has "
3980 << "cells on " << arg2 << " levels.");
3986 DeclException0(ExcGridReadError);
3991 DeclException0(ExcFacesHaveNoLevel);
3997 DeclException1(ExcEmptyLevel,
3998 int,
3999 << "You tried to do something on level " << arg1
4000 << ", but this level is empty.");
4001
4009 DeclException1(ExcBoundaryIdNotFound,
4011 << "The given boundary_id " << arg1
4012 << " is not defined in this Triangulation!");
4013
4020 ExcInconsistentCoarseningFlags,
4021 "A cell is flagged for coarsening, but either not all of its siblings "
4022 "are active or flagged for coarsening as well. Please clean up all "
4023 "coarsen flags on your triangulation via "
4024 "Triangulation::prepare_coarsening_and_refinement() beforehand!");
4025
4028protected:
4034
4039 std::vector<ReferenceCell> reference_cells;
4040
4054 static void
4055 write_bool_vector(const unsigned int magic_number1,
4056 const std::vector<bool> &v,
4057 const unsigned int magic_number2,
4058 std::ostream &out);
4059
4064 static void
4065 read_bool_vector(const unsigned int magic_number1,
4066 std::vector<bool> &v,
4067 const unsigned int magic_number2,
4068 std::istream &in);
4069
4074 void
4076
4080 virtual void
4082
4083
4084private:
4089 std::unique_ptr<
4092
4099 std::vector<GridTools::PeriodicFacePair<cell_iterator>>
4101
4106 std::map<std::pair<cell_iterator, unsigned int>,
4107 std::pair<std::pair<cell_iterator, unsigned int>, unsigned char>>
4109
4120 TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
4123 using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
4124 using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
4125 using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
4126
4137 begin_raw(const unsigned int level = 0) const;
4138
4144 end_raw(const unsigned int level) const;
4145
4160 begin_raw_line(const unsigned int level = 0) const;
4161
4180 begin_line(const unsigned int level = 0) const;
4181
4200 begin_active_line(const unsigned int level = 0) const;
4201
4207 end_line() const;
4208
4236 begin_raw_quad(const unsigned int level = 0) const;
4237
4256 begin_quad(const unsigned int level = 0) const;
4257
4276 begin_active_quad(const unsigned int level = 0) const;
4277
4283 end_quad() const;
4284
4311 begin_raw_hex(const unsigned int level = 0) const;
4312
4331 begin_hex(const unsigned int level = 0) const;
4332
4351 begin_active_hex(const unsigned int level = 0) const;
4352
4358 end_hex() const;
4359
4376 void
4378
4382 void
4384
4391 void
4393
4397 void
4399
4403 void
4405
4421
4428 void
4430
4435 void
4437
4453 virtual unsigned int
4455 const types::coarse_cell_id coarse_cell_id) const;
4456
4457
4470 virtual types::coarse_cell_id
4472 const unsigned int coarse_cell_index) const;
4473
4478 std::vector<
4479 std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
4481
4487 std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
4489
4490
4494 std::vector<Point<spacedim>> vertices;
4495
4499 std::vector<bool> vertices_used;
4500
4505 std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4507
4512
4513
4519
4530
4545 std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4547
4548
4568 std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4570
4571 // make a couple of classes friends
4572 template <int, int, int>
4573 friend class TriaAccessorBase;
4574 template <int, int, int>
4575 friend class TriaAccessor;
4576 friend class TriaAccessor<0, 1, spacedim>;
4577
4578 friend class CellAccessor<dim, spacedim>;
4579
4580 friend struct ::internal::TriaAccessorImplementation::Implementation;
4581
4582 friend struct ::internal::TriangulationImplementation::Implementation;
4583 friend struct ::internal::TriangulationImplementation::
4584 ImplementationMixedMesh;
4585
4586 friend class ::internal::TriangulationImplementation::TriaObjects;
4587
4588 // explicitly check for sensible template arguments, but not on windows
4589 // because MSVC creates bogus warnings during normal compilation
4590#ifndef DEAL_II_MSVC
4591 static_assert(dim <= spacedim,
4592 "The dimension <dim> of a Triangulation must be less than or "
4593 "equal to the space dimension <spacedim> in which it lives.");
4594#endif
4595};
4596
4597
4598#ifndef DOXYGEN
4599
4600
4601
4602namespace internal
4603{
4604 namespace TriangulationImplementation
4605 {
4606 template <class Archive>
4607 void
4608 NumberCache<1>::serialize(Archive &ar, const unsigned int)
4609 {
4610 ar &n_levels;
4611 ar &n_lines &n_lines_level;
4612 ar &n_active_lines &n_active_lines_level;
4613 }
4614
4615
4616 template <class Archive>
4617 void
4618 NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4619 {
4620 this->NumberCache<1>::serialize(ar, version);
4621
4622 ar &n_quads &n_quads_level;
4623 ar &n_active_quads &n_active_quads_level;
4624 }
4625
4626
4627 template <class Archive>
4628 void
4629 NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4630 {
4631 this->NumberCache<2>::serialize(ar, version);
4632
4633 ar &n_hexes &n_hexes_level;
4634 ar &n_active_hexes &n_active_hexes_level;
4635 }
4636
4637 } // namespace TriangulationImplementation
4638} // namespace internal
4639
4640
4641template <int dim, int spacedim>
4644 const unsigned int index) const
4645{
4646 AssertIndexRange(index, vertices_used.size());
4647 return vertices_used[index];
4648}
4649
4650
4651
4652template <int dim, int spacedim>
4654inline unsigned int Triangulation<dim, spacedim>::n_levels() const
4655{
4656 return number_cache.n_levels;
4657}
4658
4659template <int dim, int spacedim>
4661inline unsigned int Triangulation<dim, spacedim>::n_global_levels() const
4662{
4663 return number_cache.n_levels;
4664}
4665
4666
4667template <int dim, int spacedim>
4669inline unsigned int Triangulation<dim, spacedim>::n_vertices() const
4670{
4671 return vertices.size();
4672}
4673
4674
4675
4676template <int dim, int spacedim>
4678inline const std::vector<Point<spacedim>>
4680{
4681 return vertices;
4682}
4683
4684
4685template <int dim, int spacedim>
4687template <class Archive>
4688void Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4689{
4690 // as discussed in the documentation, do not store the signals as
4691 // well as boundary and manifold description but everything else
4692 ar &smooth_grid;
4693
4694 unsigned int n_levels = levels.size();
4695 ar &n_levels;
4696 for (const auto &level : levels)
4697 ar &level;
4698
4699 // boost dereferences a nullptr when serializing a nullptr
4700 // at least up to 1.65.1. This causes problems with clang-5.
4701 // Therefore, work around it.
4702 bool faces_is_nullptr = (faces.get() == nullptr);
4703 ar &faces_is_nullptr;
4704 if (!faces_is_nullptr)
4705 ar &faces;
4706
4707 ar &vertices;
4708 ar &vertices_used;
4709
4710 ar &anisotropic_refinement;
4711 ar &number_cache;
4712
4713 ar &check_for_distorted_cells;
4714
4715 if (dim == 1)
4716 {
4717 ar &vertex_to_boundary_id_map_1d;
4718 ar &vertex_to_manifold_id_map_1d;
4719 }
4720}
4721
4722
4723
4724template <int dim, int spacedim>
4726template <class Archive>
4727void Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4728{
4729 // clear previous content. this also calls the respective signal
4730 clear();
4731
4732 // as discussed in the documentation, do not store the signals as
4733 // well as boundary and manifold description but everything else
4734 ar &smooth_grid;
4735
4736 unsigned int size;
4737 ar &size;
4738 levels.resize(size);
4739 for (auto &level_ : levels)
4740 {
4741 std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4742 ar &level;
4743 level_ = std::move(level);
4744 }
4745
4746 // Workaround for nullptr, see in save().
4747 bool faces_is_nullptr = true;
4748 ar &faces_is_nullptr;
4749 if (!faces_is_nullptr)
4750 ar &faces;
4751
4752 ar &vertices;
4753 ar &vertices_used;
4754
4755 ar &anisotropic_refinement;
4756 ar &number_cache;
4757
4758 // the levels do not serialize the active_cell_indices because
4759 // they are easy enough to rebuild upon re-loading data. do
4760 // this here. don't forget to first resize the fields appropriately
4761 {
4762 for (const auto &level : levels)
4763 {
4764 level->active_cell_indices.resize(level->refine_flags.size());
4765 level->global_active_cell_indices.resize(level->refine_flags.size());
4766 level->global_level_cell_indices.resize(level->refine_flags.size());
4767 }
4768 reset_cell_vertex_indices_cache();
4769 reset_active_cell_indices();
4770 reset_global_cell_indices();
4771 }
4772
4773 reset_policy();
4774
4775 bool my_check_for_distorted_cells;
4776 ar &my_check_for_distorted_cells;
4777
4778 Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4779 ExcMessage("The triangulation loaded into here must have the "
4780 "same setting with regard to reporting distorted "
4781 "cell as the one previously stored."));
4782
4783 if (dim == 1)
4784 {
4785 ar &vertex_to_boundary_id_map_1d;
4786 ar &vertex_to_manifold_id_map_1d;
4787 }
4788
4789 // trigger the create signal to indicate
4790 // that new content has been imported into
4791 // the triangulation
4792 signals.create();
4793}
4794
4795
4796
4797template <int dim, int spacedim>
4799inline unsigned int Triangulation<dim, spacedim>::
4801 const types::coarse_cell_id coarse_cell_id) const
4802{
4803 return coarse_cell_id;
4804}
4805
4806
4807
4808template <int dim, int spacedim>
4812 const unsigned int coarse_cell_index) const
4813{
4814 return coarse_cell_index;
4815}
4816
4817
4818
4819/* -------------- declaration of explicit specializations ------------- */
4820
4821template <>
4822unsigned int
4824template <>
4825unsigned int
4826Triangulation<1, 1>::n_quads(const unsigned int level) const;
4827template <>
4828unsigned int
4829Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4830template <>
4831unsigned int
4832Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4833template <>
4834unsigned int
4835Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4836template <>
4837unsigned int
4839template <>
4840unsigned int
4841Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4842template <>
4843unsigned int
4845template <>
4846unsigned int
4847Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4848template <>
4849unsigned int
4850Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4851template <>
4852unsigned int
4854template <>
4855unsigned int
4857template <>
4858unsigned int
4859Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4860template <>
4861unsigned int
4862Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4863
4864template <>
4865unsigned int
4867
4868
4869// -------------------------------------------------------------------
4870// -- Explicit specializations for codimension one grids
4871
4872
4873template <>
4874unsigned int
4876template <>
4877unsigned int
4878Triangulation<1, 2>::n_quads(const unsigned int level) const;
4879template <>
4880unsigned int
4881Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4882template <>
4883unsigned int
4884Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4885template <>
4886unsigned int
4887Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4888template <>
4889unsigned int
4890Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4891template <>
4892unsigned int
4894template <>
4895unsigned int
4897
4898// -------------------------------------------------------------------
4899// -- Explicit specializations for codimension two grids
4900
4901template <>
4902unsigned int
4904template <>
4905unsigned int
4906Triangulation<1, 3>::n_quads(const unsigned int level) const;
4907template <>
4908unsigned int
4909Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4910template <>
4911unsigned int
4912Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4913template <>
4914unsigned int
4915Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4916template <>
4917unsigned int
4918Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4919template <>
4920unsigned int
4922template <>
4923unsigned int
4925
4926template <>
4927bool
4929template <>
4930bool
4932template <>
4933bool
4935
4936
4937extern template class Triangulation<1, 1>;
4938extern template class Triangulation<1, 2>;
4939extern template class Triangulation<1, 3>;
4940extern template class Triangulation<2, 2>;
4941extern template class Triangulation<2, 3>;
4942extern template class Triangulation<3, 3>;
4943
4944#endif // DOXYGEN
4945
4947
4948// Include tria_accessor.h here, so that it is possible for an end
4949// user to use the iterators of Triangulation<dim> directly without
4950// the need to include tria_accessor.h separately. (Otherwise the
4951// iterators are an 'opaque' or 'incomplete' type.)
4953
4954#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:4123
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:4499
bool anisotropic_refinement
Definition tria.h:4511
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, unsigned char > > periodic_face_map
Definition tria.h:4108
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:4569
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:3941
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:4494
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:4518
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:4546
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:4039
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:4488
::CellStatus CellStatus
Definition tria.h:2239
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:4124
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:4529
void load_user_flags_quad(const std::vector< bool > &v)
void save_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &file_basename) const
void load_user_flags_line(const std::vector< bool > &v)
void save_user_indices_hex(std::vector< unsigned int > &v) const
DistortedCellList execute_refinement()
void update_cell_relations()
active_line_iterator begin_active_line(const unsigned int level=0) const
void save_user_indices_quad(std::vector< unsigned int > &v) const
void load_user_pointers_hex(const std::vector< void * > &v)
void pack_data_serial()
cell_iterator end() const
virtual bool has_hanging_nodes() const
std::vector< GridTools::PeriodicFacePair< cell_iterator > > periodic_face_pairs_level_0
Definition tria.h:4100
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:3875
void load_coarsen_flags(std::istream &out)
quad_iterator end_quad() const
line_iterator begin_line(const unsigned int level=0) const
unsigned int max_adjacent_cells() const
vertex_iterator begin_vertex() const
void clear_user_flags()
unsigned int n_hexs() const
vertex_iterator end_vertex() const
bool vertex_used(const unsigned int index) const
void load_user_pointers_line(const std::vector< void * > &v)
void save_user_flags(std::vector< bool > &v) const
hex_iterator end_hex() const
hex_iterator begin_hex(const unsigned int level=0) const
virtual unsigned int n_global_levels() const
active_cell_iterator end_active(const unsigned int level) const
bool is_mixed_mesh() const
cell_iterator last() const
unsigned int n_active_quads() const
void save_coarsen_flags(std::vector< bool > &v) const
void load_user_indices_hex(const std::vector< unsigned int > &v)
unsigned int n_raw_quads() const
void save_user_pointers(std::vector< void * > &v) const
face_iterator begin_face() const
unsigned int n_cells() const
virtual bool prepare_coarsening_and_refinement()
unsigned int n_active_quads(const unsigned int level) const
void unpack_data_serial()
const std::vector< bool > & get_used_vertices() const
void load_user_flags(const std::vector< bool > &v)
typename IteratorSelector::raw_hex_iterator raw_hex_iterator
Definition tria.h:4125
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
Definition tria.h:4506
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition tria.h:4033
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition tria.h:4091
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition tria.h:2524
internal::CellAttachedDataSerializer< dim, spacedim > data_serializer
Definition tria.h:3943
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:4480
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_DEPRECATED
Definition config.h:205
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:175
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
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:466
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:534
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:489
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:511
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:2282
virtual ~DistortedCellList() noexcept override
boost::signals2::signal< void()> mesh_movement
Definition tria.h:2339
boost::signals2::signal< void()> post_distributed_load
Definition tria.h:2513
boost::signals2::signal< void()> pre_distributed_save
Definition tria.h:2491
boost::signals2::signal< unsigned int(const cell_iterator &, const ::CellStatus), CellWeightSum< unsigned int > > weight
Definition tria.h:2437
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition tria.h:2370
boost::signals2::signal< void()> pre_distributed_refinement
Definition tria.h:2450
boost::signals2::signal< void()> post_distributed_refinement
Definition tria.h:2469
boost::signals2::signal< void()> pre_refinement
Definition tria.h:2315
boost::signals2::signal< void()> any_change
Definition tria.h:2396
boost::signals2::signal< void()> create
Definition tria.h:2306
boost::signals2::signal< void()> clear
Definition tria.h:2385
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition tria.h:2360
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition tria.h:2350
boost::signals2::signal< void()> post_refinement
Definition tria.h:2322
boost::signals2::signal< void()> pre_partition
Definition tria.h:2330
boost::signals2::signal< void()> pre_distributed_repartition
Definition tria.h:2476
boost::signals2::signal< void()> pre_distributed_load
Definition tria.h:2506
boost::signals2::signal< void()> post_p4est_refinement
Definition tria.h:2459
boost::signals2::signal< void()> post_distributed_repartition
Definition tria.h:2483
boost::signals2::signal< void()> post_distributed_save
Definition tria.h:2498
std::vector< pack_callback_t > pack_callbacks_variable
Definition tria.h:362
unsigned int n_attached_deserialize
Definition tria.h:352
std::vector< pack_callback_t > pack_callbacks_fixed
Definition tria.h:361
unsigned int n_attached_data_sets
Definition tria.h:346
std::function< std::vector< char >(cell_iterator, CellStatus)> pack_callback_t
Definition tria.h:355
std::shared_ptr< const Utilities::MPI::Partitioner > active_cell_index_partitioner
Definition tria.h:184
void serialize(Archive &ar, const unsigned int version)
std::vector< unsigned int > n_active_lines_level
Definition tria.h:178
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > level_cell_index_partitioners
Definition tria.h:190
std::vector< unsigned int > n_active_quads_level
Definition tria.h:248
void serialize(Archive &ar, const unsigned int version)
void serialize(Archive &ar, const unsigned int version)
std::vector< unsigned int > n_active_hexes_level
Definition tria.h:307