Reference documentation for deal.II version GIT relicensing-35-gbe64a602e6 2024-03-01 21:20:01+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 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_tria_h
16#define dealii_tria_h
17
18
19#include <deal.II/base/config.h>
20
24#include <deal.II/base/point.h>
27
33
34#include <boost/serialization/map.hpp>
35#include <boost/serialization/split_member.hpp>
36#include <boost/serialization/unique_ptr.hpp>
37#include <boost/serialization/vector.hpp>
38#include <boost/signals2.hpp>
39
40#include <bitset>
41#include <functional>
42#include <list>
43#include <map>
44#include <memory>
45#include <numeric>
46#include <vector>
47
48
50
51#ifdef signals
52# error \
53 "The name 'signals' is already defined. You are most likely using the QT library \
54and using the 'signals' keyword. You can either #include the Qt headers (or any conflicting headers) \
55*after* the deal.II headers or you can define the 'QT_NO_KEYWORDS' macro and use the 'Q_SIGNALS' macro."
56#endif
57
58// Forward declarations
59#ifndef DOXYGEN
60template <int dim, int spacedim>
61class Manifold;
62
63template <int dim>
64struct CellData;
65
66struct SubCellData;
67
69{
70 template <int, int>
71 struct Description;
72}
73
74namespace GridTools
75{
76 template <typename CellIterator>
77 struct PeriodicFacePair;
78}
79
80template <int, int, int>
81class TriaAccessor;
82template <int spacedim>
83class TriaAccessor<0, 1, spacedim>;
84template <int, int, int>
86
87namespace internal
88{
89 namespace TriangulationImplementation
90 {
91 class TriaFaces;
92
93 class TriaObjects;
94
95 template <int, int>
96 class Policy;
97
103 struct Implementation;
104 struct ImplementationMixedMesh;
105 } // namespace TriangulationImplementation
106
107 namespace TriaAccessorImplementation
108 {
109 struct Implementation;
110 }
111} // namespace internal
112#endif
113
114
115/*------------------------------------------------------------------------*/
116
117
118namespace internal
119{
124 namespace TriangulationImplementation
125 {
137 template <int dim>
139 {};
140
152 template <>
153 struct NumberCache<1>
154 {
158 unsigned int n_levels;
159
163 unsigned int n_lines;
164
168 std::vector<unsigned int> n_lines_level;
169
173 unsigned int n_active_lines;
174
178 std::vector<unsigned int> n_active_lines_level;
179
183 std::shared_ptr<const Utilities::MPI::Partitioner>
185
189 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
191
195 NumberCache();
196
201 std::size_t
202 memory_consumption() const;
203
209 template <class Archive>
210 void
211 serialize(Archive &ar, const unsigned int version);
212 };
213
214
227 template <>
228 struct NumberCache<2> : public NumberCache<1>
229 {
233 unsigned int n_quads;
234
238 std::vector<unsigned int> n_quads_level;
239
243 unsigned int n_active_quads;
244
248 std::vector<unsigned int> n_active_quads_level;
249
253 NumberCache();
254
259 std::size_t
260 memory_consumption() const;
261
267 template <class Archive>
268 void
269 serialize(Archive &ar, const unsigned int version);
270 };
271
272
286 template <>
287 struct NumberCache<3> : public NumberCache<2>
288 {
292 unsigned int n_hexes;
293
297 std::vector<unsigned int> n_hexes_level;
298
302 unsigned int n_active_hexes;
303
307 std::vector<unsigned int> n_active_hexes_level;
308
312 NumberCache();
313
318 std::size_t
319 memory_consumption() const;
320
326 template <class Archive>
327 void
328 serialize(Archive &ar, const unsigned int version);
329 };
330 } // namespace TriangulationImplementation
331
332
336 template <int dim, int spacedim = dim>
339 {
341
347
353
355 std::function<std::vector<char>(cell_iterator, CellStatus)>;
356
361 std::vector<pack_callback_t> pack_callbacks_fixed;
362 std::vector<pack_callback_t> pack_callbacks_variable;
363 };
364
375 template <int dim, int spacedim = dim>
378 {
379 public:
381
387 using cell_relation_t = typename std::pair<cell_iterator, CellStatus>;
388
390
399 void
400 pack_data(
401 const std::vector<cell_relation_t> &cell_relations,
402 const std::vector<
404 &pack_callbacks_fixed,
405 const std::vector<
407 &pack_callbacks_variable,
408 const MPI_Comm &mpi_communicator);
409
417 void
418 unpack_cell_status(std::vector<cell_relation_t> &cell_relations) const;
419
432 void
433 unpack_data(
434 const std::vector<cell_relation_t> &cell_relations,
435 const unsigned int handle,
436 const std::function<
437 void(const cell_iterator &,
438 const CellStatus &,
439 const boost::iterator_range<std::vector<char>::const_iterator> &)>
440 &unpack_callback) const;
441
456 void
457 save(const unsigned int global_first_cell,
458 const unsigned int global_num_cells,
459 const std::string &filename,
460 const MPI_Comm &mpi_communicator) const;
461
480 void
481 load(const unsigned int global_first_cell,
482 const unsigned int global_num_cells,
483 const unsigned int local_num_cells,
484 const std::string &filename,
485 const unsigned int n_attached_deserialize_fixed,
486 const unsigned int n_attached_deserialize_variable,
487 const MPI_Comm &mpi_communicator);
488
495 void
496 clear();
497
502
513 std::vector<unsigned int> sizes_fixed_cumulative;
514
519 std::vector<char> src_data_fixed;
520 std::vector<char> dest_data_fixed;
521
526 std::vector<int> src_sizes_variable;
527 std::vector<int> dest_sizes_variable;
528 std::vector<char> src_data_variable;
529 std::vector<char> dest_data_variable;
530 };
531} // namespace internal
532
533
534/*------------------------------------------------------------------------*/
535
536
1334template <int dim, int spacedim = dim>
1337{
1338private:
1345
1346public:
1352 {
1357 none = 0x0,
1401 limit_level_difference_at_vertices = 0x1,
1422 eliminate_unrefined_islands = 0x2,
1438 patch_level_1 = 0x4,
1462 coarsest_level_1 = 0x8,
1487 allow_anisotropic_smoothing = 0x10,
1520 eliminate_refined_inner_islands = 0x100,
1525 eliminate_refined_boundary_islands = 0x200,
1531 do_not_produce_unrefined_islands = 0x400,
1532
1537 smoothing_on_refinement =
1538 (limit_level_difference_at_vertices | eliminate_unrefined_islands),
1543 smoothing_on_coarsening =
1544 (eliminate_refined_inner_islands | eliminate_refined_boundary_islands |
1545 do_not_produce_unrefined_islands),
1546
1552 maximum_smoothing = 0xffff ^ allow_anisotropic_smoothing
1554
1571
1577
1595
1609 using face_iterator = TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1610
1623 TriaActiveIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1624
1634
1649
1657 using line_iterator = typename IteratorSelector::line_iterator;
1658
1672 using active_line_iterator = typename IteratorSelector::active_line_iterator;
1673
1681 using quad_iterator = typename IteratorSelector::quad_iterator;
1682
1696 using active_quad_iterator = typename IteratorSelector::active_quad_iterator;
1697
1705 using hex_iterator = typename IteratorSelector::hex_iterator;
1706
1716 using active_hex_iterator = typename IteratorSelector::active_hex_iterator;
1717
1737 {
1744 virtual ~DistortedCellList() noexcept override;
1745
1750 std::list<typename Triangulation<dim, spacedim>::cell_iterator>
1751 distorted_cells;
1752 };
1753
1757 static constexpr unsigned int dimension = dim;
1758
1762 static constexpr unsigned int space_dimension = spacedim;
1763
1778 Triangulation(const MeshSmoothing smooth_grid = none,
1779 const bool check_for_distorted_cells = false);
1780
1796 Triangulation(const Triangulation<dim, spacedim> &) = delete;
1797
1804 Triangulation(Triangulation<dim, spacedim> &&tria) noexcept;
1805
1810 operator=(Triangulation<dim, spacedim> &&tria) noexcept;
1811
1815 virtual ~Triangulation() override;
1816
1823 virtual void
1824 clear();
1825
1830 virtual MPI_Comm
1831 get_communicator() const;
1832
1838 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1839 global_active_cell_index_partitioner() const;
1840
1846 virtual std::weak_ptr<const Utilities::MPI::Partitioner>
1847 global_level_cell_index_partitioner(const unsigned int level) const;
1848
1853 virtual void
1854 set_mesh_smoothing(const MeshSmoothing mesh_smoothing);
1855
1859 virtual const MeshSmoothing &
1860 get_mesh_smoothing() const;
1861
1884 void
1885 set_manifold(const types::manifold_id number,
1886 const Manifold<dim, spacedim> &manifold_object);
1887
1904 void
1905 reset_manifold(const types::manifold_id manifold_number);
1906
1923 void
1924 reset_all_manifolds();
1925
1934 void
1935 set_all_manifold_ids(const types::manifold_id number);
1936
1945 void
1946 set_all_manifold_ids_on_boundary(const types::manifold_id number);
1947
1957 void
1958 set_all_manifold_ids_on_boundary(const types::boundary_id b_id,
1959 const types::manifold_id number);
1960
1972 const Manifold<dim, spacedim> &
1973 get_manifold(const types::manifold_id number) const;
1974
1987 virtual std::vector<types::boundary_id>
1988 get_boundary_ids() const;
1989
2002 virtual std::vector<types::manifold_id>
2003 get_manifold_ids() const;
2004
2026 virtual void
2027 copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
2028
2077 virtual void
2078 create_triangulation(const std::vector<Point<spacedim>> &vertices,
2079 const std::vector<CellData<dim>> &cells,
2080 const SubCellData &subcelldata);
2081
2094 virtual void
2095 create_triangulation(
2096 const TriangulationDescription::Description<dim, spacedim>
2097 &construction_data);
2098
2107 void
2108 flip_all_direction_flags();
2109
2121 void
2122 set_all_refine_flags();
2123
2146 void
2147 refine_global(const unsigned int times = 1);
2148
2161 void
2162 coarsen_global(const unsigned int times = 1);
2163
2195 virtual void
2196 execute_coarsening_and_refinement();
2197
2227 virtual bool
2228 prepare_coarsening_and_refinement();
2229
2244 using CellStatus DEAL_II_DEPRECATED_EARLY = ::CellStatus;
2245
2250 static constexpr auto CELL_PERSIST DEAL_II_DEPRECATED_EARLY =
2252
2257 static constexpr auto CELL_REFINE DEAL_II_DEPRECATED_EARLY =
2259
2264 static constexpr auto CELL_COARSEN DEAL_II_DEPRECATED_EARLY =
2266
2271 static constexpr auto CELL_INVALID DEAL_II_DEPRECATED_EARLY =
2273
2274
2280 template <typename T>
2282 {
2283 using result_type = T;
2284
2285 template <typename InputIterator>
2286 T
2287 operator()(InputIterator first, InputIterator last) const
2288 {
2289 return std::accumulate(first, last, T());
2290 }
2291 };
2292
2302 struct Signals
2303 {
2311 boost::signals2::signal<void()> create;
2312
2320 boost::signals2::signal<void()> pre_refinement;
2321
2327 boost::signals2::signal<void()> post_refinement;
2328
2335 boost::signals2::signal<void()> pre_partition;
2336
2344 boost::signals2::signal<void()> mesh_movement;
2345
2353 boost::signals2::signal<void(
2354 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2356
2363 boost::signals2::signal<void(
2364 const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2366
2373 boost::signals2::signal<void(
2374 const Triangulation<dim, spacedim> &destination_tria)>
2376
2390 boost::signals2::signal<void()> clear;
2391
2401 boost::signals2::signal<void()> any_change;
2402
2439 boost::signals2::signal<unsigned int(const cell_iterator &,
2440 const ::CellStatus),
2443
2455 boost::signals2::signal<void()> pre_distributed_refinement;
2456
2464 boost::signals2::signal<void()> post_p4est_refinement;
2465
2474 boost::signals2::signal<void()> post_distributed_refinement;
2475
2481 boost::signals2::signal<void()> pre_distributed_repartition;
2482
2488 boost::signals2::signal<void()> post_distributed_repartition;
2489
2496 boost::signals2::signal<void()> pre_distributed_save;
2497
2503 boost::signals2::signal<void()> post_distributed_save;
2504
2511 boost::signals2::signal<void()> pre_distributed_load;
2512
2518 boost::signals2::signal<void()> post_distributed_load;
2519 };
2520
2530
2542 void
2543 save_refine_flags(std::ostream &out) const;
2544
2548 void
2549 save_refine_flags(std::vector<bool> &v) const;
2550
2554 void
2555 load_refine_flags(std::istream &in);
2556
2560 void
2561 load_refine_flags(const std::vector<bool> &v);
2562
2566 void
2567 save_coarsen_flags(std::ostream &out) const;
2568
2572 void
2573 save_coarsen_flags(std::vector<bool> &v) const;
2574
2578 void
2579 load_coarsen_flags(std::istream &out);
2580
2584 void
2585 load_coarsen_flags(const std::vector<bool> &v);
2586
2591 bool
2593
2605 void
2607
2613 void
2614 save_user_flags(std::ostream &out) const;
2615
2621 void
2622 save_user_flags(std::vector<bool> &v) const;
2623
2628 void
2629 load_user_flags(std::istream &in);
2630
2635 void
2636 load_user_flags(const std::vector<bool> &v);
2637
2642 void
2644
2649 void
2650 save_user_flags_line(std::ostream &out) const;
2651
2657 void
2658 save_user_flags_line(std::vector<bool> &v) const;
2659
2664 void
2665 load_user_flags_line(std::istream &in);
2666
2671 void
2672 load_user_flags_line(const std::vector<bool> &v);
2673
2678 void
2680
2685 void
2686 save_user_flags_quad(std::ostream &out) const;
2687
2693 void
2694 save_user_flags_quad(std::vector<bool> &v) const;
2695
2700 void
2701 load_user_flags_quad(std::istream &in);
2702
2707 void
2708 load_user_flags_quad(const std::vector<bool> &v);
2709
2710
2715 void
2717
2722 void
2723 save_user_flags_hex(std::ostream &out) const;
2724
2730 void
2731 save_user_flags_hex(std::vector<bool> &v) const;
2732
2737 void
2738 load_user_flags_hex(std::istream &in);
2739
2744 void
2745 load_user_flags_hex(const std::vector<bool> &v);
2746
2752 void
2754
2760 void
2761 save_user_indices(std::vector<unsigned int> &v) const;
2762
2767 void
2768 load_user_indices(const std::vector<unsigned int> &v);
2769
2775 void
2776 save_user_pointers(std::vector<void *> &v) const;
2777
2782 void
2783 load_user_pointers(const std::vector<void *> &v);
2784
2790 void
2791 save_user_indices_line(std::vector<unsigned int> &v) const;
2792
2797 void
2798 load_user_indices_line(const std::vector<unsigned int> &v);
2799
2805 void
2806 save_user_indices_quad(std::vector<unsigned int> &v) const;
2807
2812 void
2813 load_user_indices_quad(const std::vector<unsigned int> &v);
2814
2820 void
2821 save_user_indices_hex(std::vector<unsigned int> &v) const;
2822
2827 void
2828 load_user_indices_hex(const std::vector<unsigned int> &v);
2834 void
2835 save_user_pointers_line(std::vector<void *> &v) const;
2836
2841 void
2842 load_user_pointers_line(const std::vector<void *> &v);
2843
2849 void
2850 save_user_pointers_quad(std::vector<void *> &v) const;
2851
2856 void
2857 load_user_pointers_quad(const std::vector<void *> &v);
2858
2864 void
2865 save_user_pointers_hex(std::vector<void *> &v) const;
2866
2871 void
2872 load_user_pointers_hex(const std::vector<void *> &v);
2873
2899 begin(const unsigned int level = 0) const;
2900
2932 begin_active(const unsigned int level = 0) const;
2933
2939 end() const;
2940
2960 end(const unsigned int level) const;
2961
2982 end_active(const unsigned int level) const;
2983
2984
2989 last() const;
2990
2996
3005 create_cell_iterator(const CellId &cell_id) const;
3006
3017 bool
3018 contains_cell(const CellId &cell_id) const;
3038
3076
3093 cell_iterators_on_level(const unsigned int level) const;
3094
3111 active_cell_iterators_on_level(const unsigned int level) const;
3112
3115 /*-------------------------------------------------------------------------*/
3116
3126 begin_face() const;
3127
3133
3139 end_face() const;
3140
3161
3164 /*-------------------------------------------------------------------------*/
3165
3177
3185
3192 end_vertex() const;
3193
3213 unsigned int
3214 n_lines() const;
3215
3219 unsigned int
3220 n_lines(const unsigned int level) const;
3221
3225 unsigned int
3227
3231 unsigned int
3232 n_active_lines(const unsigned int level) const;
3233
3237 unsigned int
3238 n_quads() const;
3239
3243 unsigned int
3244 n_quads(const unsigned int level) const;
3245
3249 unsigned int
3251
3255 unsigned int
3256 n_active_quads(const unsigned int level) const;
3257
3261 unsigned int
3262 n_hexs() const;
3263
3268 unsigned int
3269 n_hexs(const unsigned int level) const;
3270
3274 unsigned int
3276
3281 unsigned int
3282 n_active_hexs(const unsigned int level) const;
3283
3288 unsigned int
3289 n_cells() const;
3290
3295 unsigned int
3296 n_cells(const unsigned int level) const;
3297
3302 unsigned int
3304
3314
3315
3320 unsigned int
3321 n_active_cells(const unsigned int level) const;
3322
3327 virtual types::coarse_cell_id
3329
3335 unsigned int
3336 n_faces() const;
3337
3343 unsigned int
3345
3363 unsigned int
3364 n_levels() const;
3365
3372 virtual unsigned int
3374
3384 virtual bool
3386
3394 unsigned int
3395 n_vertices() const;
3396
3405 const std::vector<Point<spacedim>> &
3407
3412 unsigned int
3414
3418 bool
3419 vertex_used(const unsigned int index) const;
3420
3425 const std::vector<bool> &
3427
3439 unsigned int
3441
3448 virtual types::subdomain_id
3450
3462
3469
3470
3487 unsigned int
3489
3499 unsigned int
3500 n_raw_lines(const unsigned int level) const;
3501
3511 unsigned int
3513
3523 unsigned int
3524 n_raw_quads(const unsigned int level) const;
3525
3535 unsigned int
3536 n_raw_hexs(const unsigned int level) const;
3537
3547 unsigned int
3548 n_raw_cells(const unsigned int level) const;
3549
3561 unsigned int
3563
3574 virtual std::size_t
3576
3586 template <class Archive>
3587 void
3588 save(Archive &ar, const unsigned int version) const;
3589
3607 template <class Archive>
3608 void
3609 load(Archive &ar, const unsigned int version);
3610
3611
3618 virtual void
3619 save(const std::string &filename) const;
3620
3624 virtual void
3625 load(const std::string &filename);
3626
3627
3643 virtual void
3645 const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3646
3650 const std::map<
3651 std::pair<cell_iterator, unsigned int>,
3652 std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3>>> &
3654
3659 const std::vector<ReferenceCell> &
3661
3666 bool
3668
3673 bool
3675
3681 bool
3683
3684#ifdef DOXYGEN
3690 template <class Archive>
3691 void
3692 serialize(Archive &archive, const unsigned int version);
3693#else
3694 // This macro defines the serialize() method that is compatible with
3695 // the templated save() and load() method that have been implemented.
3696 BOOST_SERIALIZATION_SPLIT_MEMBER()
3697#endif
3698
3703public:
3811 unsigned int
3813 const std::function<std::vector<char>(const cell_iterator &,
3814 const ::CellStatus)>
3815 &pack_callback,
3816 const bool returns_variable_size_data);
3817
3867 void
3869 const unsigned int handle,
3870 const std::function<
3871 void(const cell_iterator &,
3872 const ::CellStatus,
3873 const boost::iterator_range<std::vector<char>::const_iterator> &)>
3874 &unpack_callback);
3875
3877
3878protected:
3885 void
3886 save_attached_data(const unsigned int global_first_cell,
3887 const unsigned int global_num_cells,
3888 const std::string &filename) const;
3889
3897 void
3898 load_attached_data(const unsigned int global_first_cell,
3899 const unsigned int global_num_cells,
3900 const unsigned int local_num_cells,
3901 const std::string &filename,
3902 const unsigned int n_attached_deserialize_fixed,
3903 const unsigned int n_attached_deserialize_variable);
3904
3915 virtual void
3918
3924 std::vector<typename internal::CellAttachedDataSerializer<dim, spacedim>::
3925 cell_relation_t>
3927
3933public:
3944 DeclException2(ExcInvalidLevel,
3945 int,
3946 int,
3947 << "You are requesting information from refinement level "
3948 << arg1
3949 << " of a triangulation, but this triangulation only has "
3950 << arg2 << " refinement levels. The given level " << arg1
3951 << " must be *less* than " << arg2 << '.');
3959 ExcTriangulationNotEmpty,
3960 int,
3961 int,
3962 << "You are trying to perform an operation on a triangulation "
3963 << "that is only allowed if the triangulation is currently empty. "
3964 << "However, it currently stores " << arg1 << " vertices and has "
3965 << "cells on " << arg2 << " levels.");
3971 DeclException0(ExcGridReadError);
3976 DeclException0(ExcFacesHaveNoLevel);
3982 DeclException1(ExcEmptyLevel,
3983 int,
3984 << "You tried to do something on level " << arg1
3985 << ", but this level is empty.");
3986
3994 DeclException1(ExcBoundaryIdNotFound,
3996 << "The given boundary_id " << arg1
3997 << " is not defined in this Triangulation!");
3998
4005 ExcInconsistentCoarseningFlags,
4006 "A cell is flagged for coarsening, but either not all of its siblings "
4007 "are active or flagged for coarsening as well. Please clean up all "
4008 "coarsen flags on your triangulation via "
4009 "Triangulation::prepare_coarsening_and_refinement() beforehand!");
4010
4013protected:
4019
4024 std::vector<ReferenceCell> reference_cells;
4025
4039 static void
4040 write_bool_vector(const unsigned int magic_number1,
4041 const std::vector<bool> &v,
4042 const unsigned int magic_number2,
4043 std::ostream &out);
4044
4049 static void
4050 read_bool_vector(const unsigned int magic_number1,
4051 std::vector<bool> &v,
4052 const unsigned int magic_number2,
4053 std::istream &in);
4054
4059 void
4061
4065 virtual void
4067
4068
4069private:
4074 std::unique_ptr<
4077
4084 std::vector<GridTools::PeriodicFacePair<cell_iterator>>
4086
4091 std::map<std::pair<cell_iterator, unsigned int>,
4092 std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3>>>
4094
4105 TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
4108 using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
4109 using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
4110 using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
4111
4122 begin_raw(const unsigned int level = 0) const;
4123
4129 end_raw(const unsigned int level) const;
4130
4145 begin_raw_line(const unsigned int level = 0) const;
4146
4165 begin_line(const unsigned int level = 0) const;
4166
4185 begin_active_line(const unsigned int level = 0) const;
4186
4192 end_line() const;
4193
4221 begin_raw_quad(const unsigned int level = 0) const;
4222
4241 begin_quad(const unsigned int level = 0) const;
4242
4261 begin_active_quad(const unsigned int level = 0) const;
4262
4268 end_quad() const;
4269
4296 begin_raw_hex(const unsigned int level = 0) const;
4297
4316 begin_hex(const unsigned int level = 0) const;
4317
4336 begin_active_hex(const unsigned int level = 0) const;
4337
4343 end_hex() const;
4344
4361 void
4363
4367 void
4369
4376 void
4378
4382 void
4384
4388 void
4390
4406
4413 void
4415
4420 void
4422
4438 virtual unsigned int
4440 const types::coarse_cell_id coarse_cell_id) const;
4441
4442
4455 virtual types::coarse_cell_id
4457 const unsigned int coarse_cell_index) const;
4458
4463 std::vector<
4464 std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
4466
4472 std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
4474
4475
4479 std::vector<Point<spacedim>> vertices;
4480
4484 std::vector<bool> vertices_used;
4485
4490 std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4492
4497
4498
4504
4515
4530 std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4532
4533
4553 std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4555
4556 // make a couple of classes friends
4557 template <int, int, int>
4558 friend class TriaAccessorBase;
4559 template <int, int, int>
4560 friend class TriaAccessor;
4561 friend class TriaAccessor<0, 1, spacedim>;
4562
4563 friend class CellAccessor<dim, spacedim>;
4564
4565 friend struct ::internal::TriaAccessorImplementation::Implementation;
4566
4567 friend struct ::internal::TriangulationImplementation::Implementation;
4568 friend struct ::internal::TriangulationImplementation::
4569 ImplementationMixedMesh;
4570
4571 friend class ::internal::TriangulationImplementation::TriaObjects;
4572
4573 // explicitly check for sensible template arguments, but not on windows
4574 // because MSVC creates bogus warnings during normal compilation
4575#ifndef DEAL_II_MSVC
4576 static_assert(dim <= spacedim,
4577 "The dimension <dim> of a Triangulation must be less than or "
4578 "equal to the space dimension <spacedim> in which it lives.");
4579#endif
4580};
4581
4582
4583#ifndef DOXYGEN
4584
4585
4586
4587namespace internal
4588{
4589 namespace TriangulationImplementation
4590 {
4591 template <class Archive>
4592 void
4593 NumberCache<1>::serialize(Archive &ar, const unsigned int)
4594 {
4595 ar &n_levels;
4596 ar &n_lines &n_lines_level;
4597 ar &n_active_lines &n_active_lines_level;
4598 }
4599
4600
4601 template <class Archive>
4602 void
4603 NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4604 {
4605 this->NumberCache<1>::serialize(ar, version);
4606
4607 ar &n_quads &n_quads_level;
4608 ar &n_active_quads &n_active_quads_level;
4609 }
4610
4611
4612 template <class Archive>
4613 void
4614 NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4615 {
4616 this->NumberCache<2>::serialize(ar, version);
4617
4618 ar &n_hexes &n_hexes_level;
4619 ar &n_active_hexes &n_active_hexes_level;
4620 }
4621
4622 } // namespace TriangulationImplementation
4623} // namespace internal
4624
4625
4626template <int dim, int spacedim>
4629 const unsigned int index) const
4630{
4631 AssertIndexRange(index, vertices_used.size());
4632 return vertices_used[index];
4633}
4634
4635
4636
4637template <int dim, int spacedim>
4639inline unsigned int Triangulation<dim, spacedim>::n_levels() const
4640{
4641 return number_cache.n_levels;
4642}
4643
4644template <int dim, int spacedim>
4646inline unsigned int Triangulation<dim, spacedim>::n_global_levels() const
4647{
4648 return number_cache.n_levels;
4649}
4650
4651
4652template <int dim, int spacedim>
4654inline unsigned int Triangulation<dim, spacedim>::n_vertices() const
4655{
4656 return vertices.size();
4657}
4658
4659
4660
4661template <int dim, int spacedim>
4663inline const std::vector<Point<spacedim>>
4665{
4666 return vertices;
4667}
4668
4669
4670template <int dim, int spacedim>
4672template <class Archive>
4673void Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4674{
4675 // as discussed in the documentation, do not store the signals as
4676 // well as boundary and manifold description but everything else
4677 ar &smooth_grid;
4678
4679 unsigned int n_levels = levels.size();
4680 ar &n_levels;
4681 for (const auto &level : levels)
4682 ar &level;
4683
4684 // boost dereferences a nullptr when serializing a nullptr
4685 // at least up to 1.65.1. This causes problems with clang-5.
4686 // Therefore, work around it.
4687 bool faces_is_nullptr = (faces.get() == nullptr);
4688 ar &faces_is_nullptr;
4689 if (!faces_is_nullptr)
4690 ar &faces;
4691
4692 ar &vertices;
4693 ar &vertices_used;
4694
4695 ar &anisotropic_refinement;
4696 ar &number_cache;
4697
4698 ar &check_for_distorted_cells;
4699
4700 if (dim == 1)
4701 {
4702 ar &vertex_to_boundary_id_map_1d;
4703 ar &vertex_to_manifold_id_map_1d;
4704 }
4705}
4706
4707
4708
4709template <int dim, int spacedim>
4711template <class Archive>
4712void Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4713{
4714 // clear previous content. this also calls the respective signal
4715 clear();
4716
4717 // as discussed in the documentation, do not store the signals as
4718 // well as boundary and manifold description but everything else
4719 ar &smooth_grid;
4720
4721 unsigned int size;
4722 ar &size;
4723 levels.resize(size);
4724 for (auto &level_ : levels)
4725 {
4726 std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4727 ar &level;
4728 level_ = std::move(level);
4729 }
4730
4731 // Workaround for nullptr, see in save().
4732 bool faces_is_nullptr = true;
4733 ar &faces_is_nullptr;
4734 if (!faces_is_nullptr)
4735 ar &faces;
4736
4737 ar &vertices;
4738 ar &vertices_used;
4739
4740 ar &anisotropic_refinement;
4741 ar &number_cache;
4742
4743 // the levels do not serialize the active_cell_indices because
4744 // they are easy enough to rebuild upon re-loading data. do
4745 // this here. don't forget to first resize the fields appropriately
4746 {
4747 for (auto &level : levels)
4748 {
4749 level->active_cell_indices.resize(level->refine_flags.size());
4750 level->global_active_cell_indices.resize(level->refine_flags.size());
4751 level->global_level_cell_indices.resize(level->refine_flags.size());
4752 }
4753 reset_cell_vertex_indices_cache();
4754 reset_active_cell_indices();
4755 reset_global_cell_indices();
4756 }
4757
4758 reset_policy();
4759
4760 bool my_check_for_distorted_cells;
4761 ar &my_check_for_distorted_cells;
4762
4763 Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4764 ExcMessage("The triangulation loaded into here must have the "
4765 "same setting with regard to reporting distorted "
4766 "cell as the one previously stored."));
4767
4768 if (dim == 1)
4769 {
4770 ar &vertex_to_boundary_id_map_1d;
4771 ar &vertex_to_manifold_id_map_1d;
4772 }
4773
4774 // trigger the create signal to indicate
4775 // that new content has been imported into
4776 // the triangulation
4777 signals.create();
4778}
4779
4780
4781
4782template <int dim, int spacedim>
4784inline unsigned int Triangulation<dim, spacedim>::
4786 const types::coarse_cell_id coarse_cell_id) const
4787{
4788 return coarse_cell_id;
4789}
4790
4791
4792
4793template <int dim, int spacedim>
4797 const unsigned int coarse_cell_index) const
4798{
4799 return coarse_cell_index;
4800}
4801
4802
4803
4804/* -------------- declaration of explicit specializations ------------- */
4805
4806template <>
4807unsigned int
4809template <>
4810unsigned int
4811Triangulation<1, 1>::n_quads(const unsigned int level) const;
4812template <>
4813unsigned int
4814Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4815template <>
4816unsigned int
4817Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4818template <>
4819unsigned int
4820Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4821template <>
4822unsigned int
4824template <>
4825unsigned int
4826Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4827template <>
4828unsigned int
4830template <>
4831unsigned int
4832Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4833template <>
4834unsigned int
4835Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4836template <>
4837unsigned int
4839template <>
4840unsigned int
4842template <>
4843unsigned int
4844Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4845template <>
4846unsigned int
4847Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4848
4849template <>
4850unsigned int
4852
4853
4854// -------------------------------------------------------------------
4855// -- Explicit specializations for codimension one grids
4856
4857
4858template <>
4859unsigned int
4861template <>
4862unsigned int
4863Triangulation<1, 2>::n_quads(const unsigned int level) const;
4864template <>
4865unsigned int
4866Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4867template <>
4868unsigned int
4869Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4870template <>
4871unsigned int
4872Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4873template <>
4874unsigned int
4875Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4876template <>
4877unsigned int
4879template <>
4880unsigned int
4882
4883// -------------------------------------------------------------------
4884// -- Explicit specializations for codimension two grids
4885
4886template <>
4887unsigned int
4889template <>
4890unsigned int
4891Triangulation<1, 3>::n_quads(const unsigned int level) const;
4892template <>
4893unsigned int
4894Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4895template <>
4896unsigned int
4897Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4898template <>
4899unsigned int
4900Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4901template <>
4902unsigned int
4903Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4904template <>
4905unsigned int
4907template <>
4908unsigned int
4910
4911template <>
4912bool
4914template <>
4915bool
4917template <>
4918bool
4920
4921
4922extern template class Triangulation<1, 1>;
4923extern template class Triangulation<1, 2>;
4924extern template class Triangulation<1, 3>;
4925extern template class Triangulation<2, 2>;
4926extern template class Triangulation<2, 3>;
4927extern template class Triangulation<3, 3>;
4928
4929#endif // DOXYGEN
4930
4932
4933// Include tria_accessor.h here, so that it is possible for an end
4934// user to use the iterators of Triangulation<dim> directly without
4935// the need to include tria_accessor.h separately. (Otherwise the
4936// iterators are an 'opaque' or 'incomplete' type.)
4938
4939#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:4108
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:4484
bool anisotropic_refinement
Definition tria.h:4496
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:4554
void save_user_pointers_quad(std::vector< void * > &v) const
void save_user_flags_hex(std::ostream &out) const
void clear_user_flags_quad()
unsigned int n_faces() const
active_hex_iterator begin_active_hex(const unsigned int level=0) const
static void read_bool_vector(const unsigned int magic_number1, std::vector< bool > &v, const unsigned int magic_number2, std::istream &in)
unsigned int n_active_lines(const unsigned int level) const
bool all_reference_cells_are_hyper_cube() const
void load_user_flags_line(std::istream &in)
void clear_user_data()
raw_hex_iterator begin_raw_hex(const unsigned int level=0) const
void save_user_flags_line(std::ostream &out) const
active_cell_iterator last_active() const
void save(Archive &ar, const unsigned int version) const
std::vector< typename internal::CellAttachedDataSerializer< dim, spacedim >::cell_relation_t > local_cell_relations
Definition tria.h:3926
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:4093
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:4479
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:4503
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:4531
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:4024
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:4473
::CellStatus CellStatus
Definition tria.h:2244
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:4109
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:4514
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:4085
unsigned int n_raw_cells(const unsigned int level) const
cell_iterator end(const unsigned int level) const
bool contains_cell(const CellId &cell_id) const
virtual void save(const std::string &filename) const
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const
unsigned int n_raw_quads(const unsigned int level) const
internal::CellAttachedData< dim, spacedim > cell_attached_data
Definition tria.h:3876
void load_coarsen_flags(std::istream &out)
quad_iterator end_quad() const
line_iterator begin_line(const unsigned int level=0) const
unsigned int max_adjacent_cells() const
vertex_iterator begin_vertex() const
void save_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &filename) const
void clear_user_flags()
unsigned int n_hexs() const
vertex_iterator end_vertex() const
bool vertex_used(const unsigned int index) const
void load_user_pointers_line(const std::vector< void * > &v)
void save_user_flags(std::vector< bool > &v) const
hex_iterator end_hex() const
virtual void load(const std::string &filename)
hex_iterator begin_hex(const unsigned int level=0) const
virtual unsigned int n_global_levels() const
active_cell_iterator end_active(const unsigned int level) const
bool is_mixed_mesh() const
virtual void update_cell_relations()
Definition tria.h:3916
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:4110
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
Definition tria.h:4491
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition tria.h:4018
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition tria.h:4076
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition tria.h:2529
internal::CellAttachedDataSerializer< dim, spacedim > data_serializer
Definition tria.h:3928
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:4465
unsigned int n_raw_hexs(const unsigned int level) const
unsigned int n_lines(const unsigned int level) const
unsigned int n_active_hexs() const
virtual unsigned int coarse_cell_id_to_coarse_cell_index(const types::coarse_cell_id coarse_cell_id) const
unsigned int n_quads(const unsigned int level) const
void load_user_flags(std::istream &in)
void reset_policy()
void save_coarsen_flags(std::ostream &out) const
active_face_iterator begin_active_face() const
void clear_user_flags_line()
raw_line_iterator begin_raw_line(const unsigned int level=0) const
static void write_bool_vector(const unsigned int magic_number1, const std::vector< bool > &v, const unsigned int magic_number2, std::ostream &out)
active_cell_iterator begin_active(const unsigned int level=0) const
void execute_coarsening()
void load_attached_data(const unsigned int global_first_cell, const unsigned int global_num_cells, const unsigned int local_num_cells, const std::string &filename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable)
void save_refine_flags(std::vector< bool > &v) const
std::vector< char > dest_data_variable
Definition tria.h:529
std::vector< char > src_data_fixed
Definition tria.h:519
std::vector< char > dest_data_fixed
Definition tria.h:520
std::vector< int > src_sizes_variable
Definition tria.h:526
std::vector< unsigned int > sizes_fixed_cumulative
Definition tria.h:513
typename std::pair< cell_iterator, CellStatus > cell_relation_t
Definition tria.h:387
std::vector< int > dest_sizes_variable
Definition tria.h:527
std::vector< char > src_data_variable
Definition tria.h:528
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
Point< 3 > vertices[4]
Point< 2 > first
Definition grid_out.cc:4613
unsigned int level
Definition grid_out.cc:4616
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
IteratorRange< active_face_iterator > active_face_iterators() const
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
IteratorRange< cell_iterator > cell_iterators() const
#define DeclException0(Exception0)
Definition exceptions.h:471
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516
static ::ExceptionBase & ExcMessage(std::string arg1)
typename IteratorSelector::hex_iterator hex_iterator
Definition tria.h:1705
typename IteratorSelector::active_quad_iterator active_quad_iterator
Definition tria.h:1696
typename IteratorSelector::active_hex_iterator active_hex_iterator
Definition tria.h:1716
typename IteratorSelector::quad_iterator quad_iterator
Definition tria.h:1681
typename IteratorSelector::line_iterator line_iterator
Definition tria.h:1657
typename IteratorSelector::active_line_iterator active_line_iterator
Definition tria.h:1672
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:2287
virtual ~DistortedCellList() noexcept override
boost::signals2::signal< void()> mesh_movement
Definition tria.h:2344
boost::signals2::signal< void()> post_distributed_load
Definition tria.h:2518
boost::signals2::signal< void()> pre_distributed_save
Definition tria.h:2496
boost::signals2::signal< unsigned int(const cell_iterator &, const ::CellStatus), CellWeightSum< unsigned int > > weight
Definition tria.h:2442
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition tria.h:2375
boost::signals2::signal< void()> pre_distributed_refinement
Definition tria.h:2455
boost::signals2::signal< void()> post_distributed_refinement
Definition tria.h:2474
boost::signals2::signal< void()> pre_refinement
Definition tria.h:2320
boost::signals2::signal< void()> any_change
Definition tria.h:2401
boost::signals2::signal< void()> create
Definition tria.h:2311
boost::signals2::signal< void()> clear
Definition tria.h:2390
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition tria.h:2365
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition tria.h:2355
boost::signals2::signal< void()> post_refinement
Definition tria.h:2327
boost::signals2::signal< void()> pre_partition
Definition tria.h:2335
boost::signals2::signal< void()> pre_distributed_repartition
Definition tria.h:2481
boost::signals2::signal< void()> pre_distributed_load
Definition tria.h:2511
boost::signals2::signal< void()> post_p4est_refinement
Definition tria.h:2464
boost::signals2::signal< void()> post_distributed_repartition
Definition tria.h:2488
boost::signals2::signal< void()> post_distributed_save
Definition tria.h:2503
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