Reference documentation for deal.II version 9.3.2
\(\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\}}\)
tria.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2021 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_tria_h
17 #define dealii_tria_h
18 
19 
20 #include <deal.II/base/config.h>
21 
24 #include <deal.II/base/point.h>
27 
28 #include <deal.II/grid/cell_id.h>
32 
33 #include <boost/serialization/map.hpp>
34 #include <boost/serialization/split_member.hpp>
35 #include <boost/serialization/unique_ptr.hpp>
36 #include <boost/serialization/vector.hpp>
37 #include <boost/signals2.hpp>
38 
39 #include <bitset>
40 #include <functional>
41 #include <list>
42 #include <map>
43 #include <memory>
44 #include <numeric>
45 #include <vector>
46 
47 
49 
50 #ifdef signals
51 # error \
52  "The name 'signals' is already defined. You are most likely using the QT library \
53 and using the 'signals' keyword. You can either #include the Qt headers (or any conflicting headers) \
54 *after* the deal.II headers or you can define the 'QT_NO_KEYWORDS' macro and use the 'Q_SIGNALS' macro."
55 #endif
56 
57 // Forward declarations
58 #ifndef DOXYGEN
59 template <int dim, int spacedim>
60 class Manifold;
61 
62 template <int dim>
63 struct CellData;
64 
65 struct SubCellData;
66 
68 {
69  template <int, int>
70  struct Description;
71 }
72 
73 namespace GridTools
74 {
75  template <typename CellIterator>
76  struct PeriodicFacePair;
77 }
78 
79 template <int, int, int>
80 class TriaAccessor;
81 template <int spacedim>
82 class TriaAccessor<0, 1, spacedim>;
83 template <int, int, int>
84 class TriaAccessorBase;
85 
86 namespace internal
87 {
88  namespace TriangulationImplementation
89  {
90  class TriaFaces;
91 
92  class TriaObjects;
93 
94  template <int, int>
95  class Policy;
96 
102  struct Implementation;
103  struct ImplementationMixedMesh;
104  } // namespace TriangulationImplementation
105 
106  namespace TriaAccessorImplementation
107  {
108  struct Implementation;
109  }
110 } // namespace internal
111 
112 namespace hp
113 {
114  template <int dim, int spacedim>
115  class DoFHandler;
116 }
117 #endif
118 
119 
120 /*------------------------------------------------------------------------*/
121 
122 
123 namespace internal
124 {
129  namespace TriangulationImplementation
130  {
142  template <int dim>
143  struct NumberCache
144  {};
145 
157  template <>
158  struct NumberCache<1>
159  {
163  unsigned int n_levels;
164 
168  unsigned int n_lines;
169 
173  std::vector<unsigned int> n_lines_level;
174 
178  unsigned int n_active_lines;
179 
183  std::vector<unsigned int> n_active_lines_level;
184 
188  NumberCache();
189 
194  std::size_t
195  memory_consumption() const;
196 
202  template <class Archive>
203  void
204  serialize(Archive &ar, const unsigned int version);
205  };
206 
207 
220  template <>
221  struct NumberCache<2> : public NumberCache<1>
222  {
226  unsigned int n_quads;
227 
231  std::vector<unsigned int> n_quads_level;
232 
236  unsigned int n_active_quads;
237 
241  std::vector<unsigned int> n_active_quads_level;
242 
246  NumberCache();
247 
252  std::size_t
253  memory_consumption() const;
254 
260  template <class Archive>
261  void
262  serialize(Archive &ar, const unsigned int version);
263  };
264 
265 
279  template <>
280  struct NumberCache<3> : public NumberCache<2>
281  {
285  unsigned int n_hexes;
286 
290  std::vector<unsigned int> n_hexes_level;
291 
295  unsigned int n_active_hexes;
296 
300  std::vector<unsigned int> n_active_hexes_level;
301 
305  NumberCache();
306 
311  std::size_t
312  memory_consumption() const;
313 
319  template <class Archive>
320  void
321  serialize(Archive &ar, const unsigned int version);
322  };
323  } // namespace TriangulationImplementation
324 } // namespace internal
325 
326 
327 /*------------------------------------------------------------------------*/
328 
329 
1126 template <int dim, int spacedim = dim>
1128 {
1129 private:
1136 
1137 public:
1143  {
1148  none = 0x0,
1323 
1337 
1344  };
1345 
1362 
1368 
1386 
1400  using face_iterator = TriaIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1401 
1414  TriaActiveIterator<TriaAccessor<dim - 1, dim, spacedim>>;
1415 
1425 
1440 
1448  using line_iterator = typename IteratorSelector::line_iterator;
1449 
1463  using active_line_iterator = typename IteratorSelector::active_line_iterator;
1464 
1472  using quad_iterator = typename IteratorSelector::quad_iterator;
1473 
1487  using active_quad_iterator = typename IteratorSelector::active_quad_iterator;
1488 
1496  using hex_iterator = typename IteratorSelector::hex_iterator;
1497 
1507  using active_hex_iterator = typename IteratorSelector::active_hex_iterator;
1508 
1528  {
1535  virtual ~DistortedCellList() noexcept override;
1536 
1541  std::list<typename Triangulation<dim, spacedim>::cell_iterator>
1543  };
1544 
1548  static const unsigned int dimension = dim;
1549 
1553  static const unsigned int space_dimension = spacedim;
1554 
1570  const bool check_for_distorted_cells = false);
1571 
1587  Triangulation(const Triangulation<dim, spacedim> &) = delete;
1588 
1595  Triangulation(Triangulation<dim, spacedim> &&tria) noexcept;
1596 
1600  Triangulation &
1601  operator=(Triangulation<dim, spacedim> &&tria) noexcept;
1602 
1606  virtual ~Triangulation() override;
1607 
1614  virtual void
1616 
1621  virtual MPI_Comm
1623 
1629  virtual void
1630  set_mesh_smoothing(const MeshSmoothing mesh_smoothing);
1631 
1635  virtual const MeshSmoothing &
1637 
1660  void
1662  const Manifold<dim, spacedim> &manifold_object);
1663 
1676  void
1677  reset_manifold(const types::manifold_id manifold_number);
1678 
1690  void
1692 
1701  void
1703 
1712  void
1714 
1724  void
1726  const types::manifold_id number);
1727 
1739  const Manifold<dim, spacedim> &
1740  get_manifold(const types::manifold_id number) const;
1741 
1754  virtual std::vector<types::boundary_id>
1756 
1769  virtual std::vector<types::manifold_id>
1771 
1798  virtual void
1799  copy_triangulation(const Triangulation<dim, spacedim> &other_tria);
1800 
1848  virtual void
1849  create_triangulation(const std::vector<Point<spacedim>> &vertices,
1850  const std::vector<CellData<dim>> & cells,
1851  const SubCellData & subcelldata);
1852 
1865  virtual void
1867  const TriangulationDescription::Description<dim, spacedim>
1868  &construction_data);
1869 
1879  virtual void
1881  const std::vector<Point<spacedim>> &vertices,
1882  const std::vector<CellData<dim>> & cells,
1883  const SubCellData & subcelldata);
1884 
1891  void
1893 
1905  void
1907 
1930  void
1931  refine_global(const unsigned int times = 1);
1932 
1945  void
1946  coarsen_global(const unsigned int times = 1);
1947 
1979  virtual void
1981 
2011  virtual bool
2013 
2014  /*
2015  * @}
2016  */
2017 
2032  {
2049  CELL_INVALID
2050  };
2051 
2057  template <typename T>
2059  {
2060  using result_type = T;
2061 
2062  template <typename InputIterator>
2063  T
2064  operator()(InputIterator first, InputIterator last) const
2065  {
2066  return std::accumulate(first, last, T());
2067  }
2068  };
2069 
2079  struct Signals
2080  {
2088  boost::signals2::signal<void()> create;
2089 
2097  boost::signals2::signal<void()> pre_refinement;
2098 
2104  boost::signals2::signal<void()> post_refinement;
2105 
2112  boost::signals2::signal<void()> pre_partition;
2113 
2121  boost::signals2::signal<void()> mesh_movement;
2122 
2130  boost::signals2::signal<void(
2131  const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2133 
2140  boost::signals2::signal<void(
2141  const typename Triangulation<dim, spacedim>::cell_iterator &cell)>
2143 
2150  boost::signals2::signal<void(
2151  const Triangulation<dim, spacedim> &destination_tria)>
2153 
2167  boost::signals2::signal<void()> clear;
2168 
2178  boost::signals2::signal<void()> any_change;
2179 
2205  boost::signals2::signal<unsigned int(const cell_iterator &,
2206  const CellStatus),
2209 
2221  boost::signals2::signal<void()> pre_distributed_refinement;
2222 
2230  boost::signals2::signal<void()> post_p4est_refinement;
2231 
2240  boost::signals2::signal<void()> post_distributed_refinement;
2241 
2252  boost::signals2::signal<void()> pre_distributed_repartition;
2253 
2259  boost::signals2::signal<void()> post_distributed_repartition;
2260 
2267  boost::signals2::signal<void()> pre_distributed_save;
2268 
2274  boost::signals2::signal<void()> post_distributed_save;
2275 
2282  boost::signals2::signal<void()> pre_distributed_load;
2283 
2289  boost::signals2::signal<void()> post_distributed_load;
2290  };
2291 
2295  mutable Signals signals;
2296 
2297  /*
2298  * @}
2299  */
2300 
2310  void
2311  save_refine_flags(std::ostream &out) const;
2312 
2316  void
2317  save_refine_flags(std::vector<bool> &v) const;
2318 
2322  void
2323  load_refine_flags(std::istream &in);
2324 
2328  void
2329  load_refine_flags(const std::vector<bool> &v);
2330 
2334  void
2335  save_coarsen_flags(std::ostream &out) const;
2336 
2340  void
2341  save_coarsen_flags(std::vector<bool> &v) const;
2342 
2346  void
2347  load_coarsen_flags(std::istream &out);
2348 
2352  void
2353  load_coarsen_flags(const std::vector<bool> &v);
2354 
2359  bool
2361 
2362  /*
2363  * @}
2364  */
2365 
2375  void
2377 
2383  void
2384  save_user_flags(std::ostream &out) const;
2385 
2391  void
2392  save_user_flags(std::vector<bool> &v) const;
2393 
2398  void
2399  load_user_flags(std::istream &in);
2400 
2405  void
2406  load_user_flags(const std::vector<bool> &v);
2407 
2412  void
2414 
2419  void
2420  save_user_flags_line(std::ostream &out) const;
2421 
2427  void
2428  save_user_flags_line(std::vector<bool> &v) const;
2429 
2434  void
2435  load_user_flags_line(std::istream &in);
2436 
2441  void
2442  load_user_flags_line(const std::vector<bool> &v);
2443 
2448  void
2450 
2455  void
2456  save_user_flags_quad(std::ostream &out) const;
2457 
2463  void
2464  save_user_flags_quad(std::vector<bool> &v) const;
2465 
2470  void
2471  load_user_flags_quad(std::istream &in);
2472 
2477  void
2478  load_user_flags_quad(const std::vector<bool> &v);
2479 
2480 
2485  void
2487 
2492  void
2493  save_user_flags_hex(std::ostream &out) const;
2494 
2500  void
2501  save_user_flags_hex(std::vector<bool> &v) const;
2502 
2507  void
2508  load_user_flags_hex(std::istream &in);
2509 
2514  void
2515  load_user_flags_hex(const std::vector<bool> &v);
2516 
2522  void
2524 
2530  void
2531  save_user_indices(std::vector<unsigned int> &v) const;
2532 
2537  void
2538  load_user_indices(const std::vector<unsigned int> &v);
2539 
2545  void
2546  save_user_pointers(std::vector<void *> &v) const;
2547 
2552  void
2553  load_user_pointers(const std::vector<void *> &v);
2554 
2560  void
2561  save_user_indices_line(std::vector<unsigned int> &v) const;
2562 
2567  void
2568  load_user_indices_line(const std::vector<unsigned int> &v);
2569 
2575  void
2576  save_user_indices_quad(std::vector<unsigned int> &v) const;
2577 
2582  void
2583  load_user_indices_quad(const std::vector<unsigned int> &v);
2584 
2590  void
2591  save_user_indices_hex(std::vector<unsigned int> &v) const;
2592 
2597  void
2598  load_user_indices_hex(const std::vector<unsigned int> &v);
2604  void
2605  save_user_pointers_line(std::vector<void *> &v) const;
2606 
2611  void
2612  load_user_pointers_line(const std::vector<void *> &v);
2613 
2619  void
2620  save_user_pointers_quad(std::vector<void *> &v) const;
2621 
2626  void
2627  load_user_pointers_quad(const std::vector<void *> &v);
2628 
2634  void
2635  save_user_pointers_hex(std::vector<void *> &v) const;
2636 
2641  void
2642  load_user_pointers_hex(const std::vector<void *> &v);
2643 
2644  /*
2645  * @}
2646  */
2647 
2671  begin(const unsigned int level = 0) const;
2672 
2704  begin_active(const unsigned int level = 0) const;
2705 
2711  end() const;
2712 
2732  end(const unsigned int level) const;
2733 
2754  end_active(const unsigned int level) const;
2755 
2756 
2761  last() const;
2762 
2767  last_active() const;
2768 
2782  create_cell_iterator(const CellId &cell_id) const;
2783 
2800 
2838 
2855  cell_iterators_on_level(const unsigned int level) const;
2856 
2873  active_cell_iterators_on_level(const unsigned int level) const;
2874 
2875  /*
2876  * @}
2877  */
2878 
2879  /*-------------------------------------------------------------------------*/
2880 
2890  begin_face() const;
2891 
2897 
2903  end_face() const;
2904 
2925 
2926  /*
2927  * @}
2928  */
2929 
2930  /*-------------------------------------------------------------------------*/
2931 
2942  begin_vertex() const;
2943 
2951 
2958  end_vertex() const;
2959 
2960  /*
2961  * @}
2962  */
2963 
2981  unsigned int
2982  n_lines() const;
2983 
2987  unsigned int
2988  n_lines(const unsigned int level) const;
2989 
2993  unsigned int
2995 
2999  unsigned int
3000  n_active_lines(const unsigned int level) const;
3001 
3005  unsigned int
3006  n_quads() const;
3007 
3011  unsigned int
3012  n_quads(const unsigned int level) const;
3013 
3017  unsigned int
3019 
3023  unsigned int
3024  n_active_quads(const unsigned int level) const;
3025 
3029  unsigned int
3030  n_hexs() const;
3031 
3036  unsigned int
3037  n_hexs(const unsigned int level) const;
3038 
3042  unsigned int
3043  n_active_hexs() const;
3044 
3049  unsigned int
3050  n_active_hexs(const unsigned int level) const;
3051 
3056  unsigned int
3057  n_cells() const;
3058 
3063  unsigned int
3064  n_cells(const unsigned int level) const;
3065 
3070  unsigned int
3072 
3080  virtual types::global_cell_index
3082 
3083 
3088  unsigned int
3089  n_active_cells(const unsigned int level) const;
3090 
3096  unsigned int
3097  n_faces() const;
3098 
3104  unsigned int
3106 
3124  unsigned int
3125  n_levels() const;
3126 
3133  virtual unsigned int
3135 
3145  virtual bool
3147 
3155  unsigned int
3156  n_vertices() const;
3157 
3166  const std::vector<Point<spacedim>> &
3167  get_vertices() const;
3168 
3173  unsigned int
3175 
3179  bool
3180  vertex_used(const unsigned int index) const;
3181 
3186  const std::vector<bool> &
3188 
3200  unsigned int
3202 
3209  virtual types::subdomain_id
3211 
3223 
3230 
3231 
3232  /*
3233  * @}
3234  */
3235 
3250  unsigned int
3251  n_raw_lines() const;
3252 
3262  unsigned int
3263  n_raw_lines(const unsigned int level) const;
3264 
3274  unsigned int
3275  n_raw_quads() const;
3276 
3286  unsigned int
3287  n_raw_quads(const unsigned int level) const;
3288 
3298  unsigned int
3299  n_raw_hexs(const unsigned int level) const;
3300 
3310  unsigned int
3311  n_raw_cells(const unsigned int level) const;
3312 
3324  unsigned int
3325  n_raw_faces() const;
3326 
3327  /*
3328  * @}
3329  */
3330 
3339  virtual std::size_t
3341 
3351  template <class Archive>
3352  void
3353  save(Archive &ar, const unsigned int version) const;
3354 
3372  template <class Archive>
3373  void
3374  load(Archive &ar, const unsigned int version);
3375 
3376 
3392  virtual void
3394  const std::vector<GridTools::PeriodicFacePair<cell_iterator>> &);
3395 
3399  const std::map<
3400  std::pair<cell_iterator, unsigned int>,
3401  std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3>>> &
3403 
3408  const std::vector<ReferenceCell> &
3410 
3415  bool
3417 
3418 #ifdef DOXYGEN
3424  template <class Archive>
3425  void
3426  serialize(Archive &archive, const unsigned int version);
3427 #else
3428  // This macro defines the serialize() method that is compatible with
3429  // the templated save() and load() method that have been implemented.
3430  BOOST_SERIALIZATION_SPLIT_MEMBER()
3431 #endif
3432 
3444  int,
3445  int,
3446  << "You are requesting information from refinement level "
3447  << arg1
3448  << " of a triangulation, but this triangulation only has "
3449  << arg2 << " refinement levels. The given level " << arg1
3450  << " must be *less* than " << arg2 << ".");
3459  int,
3460  int,
3461  << "You are trying to perform an operation on a triangulation "
3462  << "that is only allowed if the triangulation is currently empty. "
3463  << "However, it currently stores " << arg1 << " vertices and has "
3464  << "cells on " << arg2 << " levels.");
3482  int,
3483  << "You tried to do something on level " << arg1
3484  << ", but this level is empty.");
3491 
3501  << "The given boundary_id " << arg1
3502  << " is not defined in this Triangulation!");
3503 
3511  "A cell is flagged for coarsening, but either not all of its siblings "
3512  "are active or flagged for coarsening as well. Please clean up all "
3513  "coarsen flags on your triangulation via "
3514  "Triangulation::prepare_coarsening_and_refinement() beforehand!");
3515 
3516  /*
3517  * @}
3518  */
3519 
3520 protected:
3526 
3531  std::vector<ReferenceCell> reference_cells;
3532 
3546  static void
3547  write_bool_vector(const unsigned int magic_number1,
3548  const std::vector<bool> &v,
3549  const unsigned int magic_number2,
3550  std::ostream & out);
3551 
3556  static void
3557  read_bool_vector(const unsigned int magic_number1,
3558  std::vector<bool> &v,
3559  const unsigned int magic_number2,
3560  std::istream & in);
3561 
3566  void
3568 
3572  virtual void
3574 
3575 
3576 private:
3581  std::unique_ptr<
3584 
3591  std::vector<GridTools::PeriodicFacePair<cell_iterator>>
3593 
3598  std::map<std::pair<cell_iterator, unsigned int>,
3599  std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3>>>
3601 
3617  TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim>>;
3620  using raw_line_iterator = typename IteratorSelector::raw_line_iterator;
3621  using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator;
3622  using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator;
3623 
3629  begin_raw(const unsigned int level = 0) const;
3630 
3636  end_raw(const unsigned int level) const;
3637 
3638  /*
3639  * @}
3640  */
3641 
3654  begin_raw_line(const unsigned int level = 0) const;
3655 
3674  begin_line(const unsigned int level = 0) const;
3675 
3694  begin_active_line(const unsigned int level = 0) const;
3695 
3701  end_line() const;
3702 
3703  /*
3704  * @}
3705  */
3706 
3732  begin_raw_quad(const unsigned int level = 0) const;
3733 
3752  begin_quad(const unsigned int level = 0) const;
3753 
3772  begin_active_quad(const unsigned int level = 0) const;
3773 
3779  end_quad() const;
3780 
3781  /*
3782  * @}
3783  */
3784 
3809  begin_raw_hex(const unsigned int level = 0) const;
3810 
3828  hex_iterator
3829  begin_hex(const unsigned int level = 0) const;
3830 
3849  begin_active_hex(const unsigned int level = 0) const;
3850 
3855  hex_iterator
3856  end_hex() const;
3857 
3858  /*
3859  * @}
3860  */
3861 
3862 
3876  void
3878 
3882  void
3884 
3891  void
3893 
3897  void
3899 
3903  void
3905 
3921 
3928  void
3930 
3935  void
3937 
3951  virtual unsigned int
3954 
3955 
3968  virtual types::coarse_cell_id
3970  const unsigned int coarse_cell_index) const;
3971 
3976  std::vector<
3977  std::unique_ptr<::internal::TriangulationImplementation::TriaLevel>>
3979 
3985  std::unique_ptr<::internal::TriangulationImplementation::TriaFaces>
3987 
3988 
3992  std::vector<Point<spacedim>> vertices;
3993 
3997  std::vector<bool> vertices_used;
3998 
4003  std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim>>>
4005 
4010 
4011 
4017 
4028 
4043  std::unique_ptr<std::map<unsigned int, types::boundary_id>>
4045 
4046 
4066  std::unique_ptr<std::map<unsigned int, types::manifold_id>>
4068 
4069  // make a couple of classes friends
4070  template <int, int, int>
4071  friend class TriaAccessorBase;
4072  template <int, int, int>
4073  friend class TriaAccessor;
4074  friend class TriaAccessor<0, 1, spacedim>;
4075 
4076  friend class CellAccessor<dim, spacedim>;
4077 
4078  friend struct ::internal::TriaAccessorImplementation::Implementation;
4079 
4080  friend struct ::internal::TriangulationImplementation::Implementation;
4081  friend struct ::internal::TriangulationImplementation::
4082  ImplementationMixedMesh;
4083 
4084  friend class ::internal::TriangulationImplementation::TriaObjects;
4085 
4086  // explicitly check for sensible template arguments, but not on windows
4087  // because MSVC creates bogus warnings during normal compilation
4088 #ifndef DEAL_II_MSVC
4089  static_assert(dim <= spacedim,
4090  "The dimension <dim> of a Triangulation must be less than or "
4091  "equal to the space dimension <spacedim> in which it lives.");
4092 #endif
4093 };
4094 
4095 
4096 #ifndef DOXYGEN
4097 
4098 
4099 
4100 namespace internal
4101 {
4102  namespace TriangulationImplementation
4103  {
4104  template <class Archive>
4105  void
4106  NumberCache<1>::serialize(Archive &ar, const unsigned int)
4107  {
4108  ar &n_levels;
4109  ar &n_lines &n_lines_level;
4110  ar &n_active_lines &n_active_lines_level;
4111  }
4112 
4113 
4114  template <class Archive>
4115  void
4116  NumberCache<2>::serialize(Archive &ar, const unsigned int version)
4117  {
4118  this->NumberCache<1>::serialize(ar, version);
4119 
4120  ar &n_quads &n_quads_level;
4121  ar &n_active_quads &n_active_quads_level;
4122  }
4123 
4124 
4125  template <class Archive>
4126  void
4127  NumberCache<3>::serialize(Archive &ar, const unsigned int version)
4128  {
4129  this->NumberCache<2>::serialize(ar, version);
4130 
4131  ar &n_hexes &n_hexes_level;
4132  ar &n_active_hexes &n_active_hexes_level;
4133  }
4134 
4135  } // namespace TriangulationImplementation
4136 } // namespace internal
4137 
4138 
4139 template <int dim, int spacedim>
4140 inline bool
4141 Triangulation<dim, spacedim>::vertex_used(const unsigned int index) const
4142 {
4143  AssertIndexRange(index, vertices_used.size());
4144  return vertices_used[index];
4145 }
4146 
4147 
4148 
4149 template <int dim, int spacedim>
4150 inline unsigned int
4152 {
4153  return number_cache.n_levels;
4154 }
4155 
4156 template <int dim, int spacedim>
4157 inline unsigned int
4159 {
4160  return number_cache.n_levels;
4161 }
4162 
4163 
4164 template <int dim, int spacedim>
4165 inline unsigned int
4167 {
4168  return vertices.size();
4169 }
4170 
4171 
4172 
4173 template <int dim, int spacedim>
4174 inline const std::vector<Point<spacedim>> &
4176 {
4177  return vertices;
4178 }
4179 
4180 
4181 template <int dim, int spacedim>
4182 template <class Archive>
4183 void
4184 Triangulation<dim, spacedim>::save(Archive &ar, const unsigned int) const
4185 {
4186  // as discussed in the documentation, do not store the signals as
4187  // well as boundary and manifold description but everything else
4188  ar &smooth_grid;
4189 
4190  unsigned int n_levels = levels.size();
4191  ar & n_levels;
4192  for (const auto &level : levels)
4193  ar &level;
4194 
4195  // boost dereferences a nullptr when serializing a nullptr
4196  // at least up to 1.65.1. This causes problems with clang-5.
4197  // Therefore, work around it.
4198  bool faces_is_nullptr = (faces.get() == nullptr);
4199  ar & faces_is_nullptr;
4200  if (!faces_is_nullptr)
4201  ar &faces;
4202 
4203  ar &vertices;
4204  ar &vertices_used;
4205 
4207  ar &number_cache;
4208 
4210 
4211  if (dim == 1)
4212  {
4215  }
4216 }
4217 
4218 
4219 
4220 template <int dim, int spacedim>
4221 template <class Archive>
4222 void
4223 Triangulation<dim, spacedim>::load(Archive &ar, const unsigned int)
4224 {
4225  // clear previous content. this also calls the respective signal
4226  clear();
4227 
4228  // as discussed in the documentation, do not store the signals as
4229  // well as boundary and manifold description but everything else
4230  ar &smooth_grid;
4231 
4232  unsigned int size;
4233  ar & size;
4234  levels.resize(size);
4235  for (auto &level_ : levels)
4236  {
4237  std::unique_ptr<internal::TriangulationImplementation::TriaLevel> level;
4238  ar & level;
4239  level_ = std::move(level);
4240  }
4241 
4242  // Workaround for nullptr, see in save().
4243  bool faces_is_nullptr = true;
4244  ar & faces_is_nullptr;
4245  if (!faces_is_nullptr)
4246  ar &faces;
4247 
4248  ar &vertices;
4249  ar &vertices_used;
4250 
4252  ar &number_cache;
4253 
4254  // the levels do not serialize the active_cell_indices because
4255  // they are easy enough to rebuild upon re-loading data. do
4256  // this here. don't forget to first resize the fields appropriately
4257  {
4258  for (auto &level : levels)
4259  {
4260  level->active_cell_indices.resize(level->refine_flags.size());
4261  level->global_active_cell_indices.resize(level->refine_flags.size());
4262  level->global_level_cell_indices.resize(level->refine_flags.size());
4263  }
4267  }
4268 
4269  reset_policy();
4270 
4271  bool my_check_for_distorted_cells;
4272  ar & my_check_for_distorted_cells;
4273 
4274  Assert(my_check_for_distorted_cells == check_for_distorted_cells,
4275  ExcMessage("The triangulation loaded into here must have the "
4276  "same setting with regard to reporting distorted "
4277  "cell as the one previously stored."));
4278 
4279  if (dim == 1)
4280  {
4283  }
4284 
4285  // trigger the create signal to indicate
4286  // that new content has been imported into
4287  // the triangulation
4288  signals.create();
4289 }
4290 
4291 
4292 
4293 template <int dim, int spacedim>
4294 inline unsigned int
4297 {
4298  return coarse_cell_id;
4299 }
4300 
4301 
4302 
4303 template <int dim, int spacedim>
4304 inline types::coarse_cell_id
4306  const unsigned int coarse_cell_index) const
4307 {
4308  return coarse_cell_index;
4309 }
4310 
4311 
4312 
4313 /* -------------- declaration of explicit specializations ------------- */
4314 
4315 template <>
4316 unsigned int
4318 template <>
4319 unsigned int
4320 Triangulation<1, 1>::n_quads(const unsigned int level) const;
4321 template <>
4322 unsigned int
4323 Triangulation<1, 1>::n_raw_quads(const unsigned int level) const;
4324 template <>
4325 unsigned int
4326 Triangulation<2, 2>::n_raw_quads(const unsigned int level) const;
4327 template <>
4328 unsigned int
4329 Triangulation<3, 3>::n_raw_quads(const unsigned int level) const;
4330 template <>
4331 unsigned int
4333 template <>
4334 unsigned int
4335 Triangulation<1, 1>::n_active_quads(const unsigned int level) const;
4336 template <>
4337 unsigned int
4339 template <>
4340 unsigned int
4341 Triangulation<1, 1>::n_raw_hexs(const unsigned int level) const;
4342 template <>
4343 unsigned int
4344 Triangulation<3, 3>::n_raw_hexs(const unsigned int level) const;
4345 template <>
4346 unsigned int
4348 template <>
4349 unsigned int
4351 template <>
4352 unsigned int
4353 Triangulation<3, 3>::n_active_hexs(const unsigned int) const;
4354 template <>
4355 unsigned int
4356 Triangulation<3, 3>::n_hexs(const unsigned int level) const;
4357 
4358 template <>
4359 unsigned int
4361 
4362 
4363 // -------------------------------------------------------------------
4364 // -- Explicit specializations for codimension one grids
4365 
4366 
4367 template <>
4368 unsigned int
4370 template <>
4371 unsigned int
4372 Triangulation<1, 2>::n_quads(const unsigned int level) const;
4373 template <>
4374 unsigned int
4375 Triangulation<1, 2>::n_raw_quads(const unsigned int level) const;
4376 template <>
4377 unsigned int
4378 Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4379 template <>
4380 unsigned int
4381 Triangulation<1, 2>::n_raw_hexs(const unsigned int level) const;
4382 template <>
4383 unsigned int
4384 Triangulation<1, 2>::n_active_quads(const unsigned int level) const;
4385 template <>
4386 unsigned int
4388 template <>
4389 unsigned int
4391 
4392 // -------------------------------------------------------------------
4393 // -- Explicit specializations for codimension two grids
4394 
4395 template <>
4396 unsigned int
4398 template <>
4399 unsigned int
4400 Triangulation<1, 3>::n_quads(const unsigned int level) const;
4401 template <>
4402 unsigned int
4403 Triangulation<1, 3>::n_raw_quads(const unsigned int level) const;
4404 template <>
4405 unsigned int
4406 Triangulation<2, 3>::n_raw_quads(const unsigned int level) const;
4407 template <>
4408 unsigned int
4409 Triangulation<1, 3>::n_raw_hexs(const unsigned int level) const;
4410 template <>
4411 unsigned int
4412 Triangulation<1, 3>::n_active_quads(const unsigned int level) const;
4413 template <>
4414 unsigned int
4416 template <>
4417 unsigned int
4419 
4420 template <>
4421 bool
4423 template <>
4424 bool
4426 template <>
4427 bool
4429 
4430 
4431 extern template class Triangulation<1, 1>;
4432 extern template class Triangulation<1, 2>;
4433 extern template class Triangulation<1, 3>;
4434 extern template class Triangulation<2, 2>;
4435 extern template class Triangulation<2, 3>;
4436 extern template class Triangulation<3, 3>;
4437 
4438 #endif // DOXYGEN
4439 
4441 
4442 // Include tria_accessor.h here, so that it is possible for an end
4443 // user to use the iterators of Triangulation<dim> directly without
4444 // the need to include tria_accessor.h separately. (Otherwise the
4445 // iterators are an 'opaque' or 'incomplete' type.)
4447 
4448 #endif
Definition: cell_id.h:71
Definition: point.h:111
unsigned int n_hexs(const unsigned int level) const
void save_user_flags_line(std::vector< bool > &v) const
virtual types::global_cell_index n_global_active_cells() const
virtual std::vector< types::boundary_id > get_boundary_ids() const
@ patch_level_1
Definition: tria.h:1229
@ smoothing_on_refinement
Definition: tria.h:1328
@ eliminate_refined_boundary_islands
Definition: tria.h:1316
@ do_not_produce_unrefined_islands
Definition: tria.h:1322
@ smoothing_on_coarsening
Definition: tria.h:1334
@ limit_level_difference_at_vertices
Definition: tria.h:1192
@ allow_anisotropic_smoothing
Definition: tria.h:1278
@ coarsest_level_1
Definition: tria.h:1253
@ eliminate_unrefined_islands
Definition: tria.h:1213
@ eliminate_refined_inner_islands
Definition: tria.h:1311
@ maximum_smoothing
Definition: tria.h:1343
virtual const MeshSmoothing & get_mesh_smoothing() const
const std::vector< bool > & get_used_vertices() const
quad_iterator begin_quad(const unsigned int level=0) const
typename IteratorSelector::raw_line_iterator raw_line_iterator
Definition: tria.h:3620
active_vertex_iterator begin_active_vertex() const
virtual MPI_Comm get_communicator() 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:3997
virtual void clear()
bool anisotropic_refinement
Definition: tria.h:4009
active_quad_iterator begin_active_quad(const unsigned int level=0) const
bool get_anisotropic_refinement_flag() const
virtual void create_triangulation(const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata)
virtual void copy_triangulation(const Triangulation< dim, spacedim > &other_tria)
std::unique_ptr< std::map< unsigned int, types::manifold_id > > vertex_to_manifold_id_map_1d
Definition: tria.h:4067
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
@ CELL_COARSEN
Definition: tria.h:2045
@ CELL_REFINE
Definition: tria.h:2041
@ CELL_INVALID
Definition: tria.h:2049
@ CELL_PERSIST
Definition: tria.h:2037
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
static const unsigned int space_dimension
Definition: tria.h:1553
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:3600
void fix_coarsen_flags()
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)
virtual void set_mesh_smoothing(const MeshSmoothing mesh_smoothing)
unsigned int n_raw_lines() const
virtual std::size_t memory_consumption() const
std::vector< Point< spacedim > > vertices
Definition: tria.h:3992
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:4016
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:4044
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:3531
void update_periodic_face_map()
unsigned int n_raw_lines(const unsigned int level) const
void clear_despite_subscriptions()
void coarsen_global(const unsigned int times=1)
void save_user_flags(std::ostream &out) const
void refine_global(const unsigned int times=1)
void load_user_flags_hex(std::istream &in)
void load_user_pointers_quad(const std::vector< void * > &v)
std::unique_ptr<::internal::TriangulationImplementation::TriaFaces > faces
Definition: tria.h:3986
virtual void add_periodicity(const std::vector< GridTools::PeriodicFacePair< cell_iterator >> &)
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
typename IteratorSelector::raw_quad_iterator raw_quad_iterator
Definition: tria.h:3621
void load_user_pointers(const std::vector< void * > &v)
::internal::TriangulationImplementation::NumberCache< dim > number_cache
Definition: tria.h:4027
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:3592
unsigned int n_raw_cells(const unsigned int level) const
cell_iterator end(const unsigned int level) const
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id(const unsigned int coarse_cell_index) const
unsigned int n_raw_quads(const unsigned int level) const
void load_coarsen_flags(std::istream &out)
quad_iterator end_quad() const
line_iterator begin_line(const unsigned int level=0) const
unsigned int max_adjacent_cells() const
vertex_iterator begin_vertex() const
void clear_user_flags()
unsigned int n_hexs() const
vertex_iterator end_vertex() const
bool vertex_used(const unsigned int index) const
void load_user_pointers_line(const std::vector< void * > &v)
void save_user_flags(std::vector< bool > &v) const
hex_iterator end_hex() const
const Triangulation< dim, spacedim > & get_triangulation() const
hex_iterator begin_hex(const unsigned int level=0) const
virtual void execute_coarsening_and_refinement()
virtual unsigned int n_global_levels() const
active_cell_iterator end_active(const unsigned int level) 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)
virtual void create_triangulation_compatibility(const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata)
unsigned int n_raw_quads() const
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifold
Definition: tria.h:4004
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 load_user_flags(const std::vector< bool > &v)
typename IteratorSelector::raw_hex_iterator raw_hex_iterator
Definition: tria.h:3622
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
void serialize(Archive &archive, const unsigned int version)
void load_refine_flags(const std::vector< bool > &v)
MeshSmoothing smooth_grid
Definition: tria.h:3525
void save_refine_flags(std::ostream &out) const
std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
Definition: tria.h:3583
Triangulation< dim, spacedim > & get_triangulation()
void save_user_flags_quad(std::ostream &out) const
Signals signals
Definition: tria.h:2295
unsigned int n_active_hexs(const unsigned int level) const
unsigned int n_vertices() const
void load(Archive &ar, const unsigned int version)
void save_user_indices(std::vector< unsigned int > &v) const
void save_user_flags_hex(std::vector< bool > &v) const
std::vector< std::unique_ptr<::internal::TriangulationImplementation::TriaLevel > > levels
Definition: tria.h:3978
unsigned int n_raw_hexs(const unsigned int level) const
void set_all_refine_flags()
unsigned int n_lines(const unsigned int level) const
unsigned int n_active_hexs() const
static const unsigned int dimension
Definition: tria.h:1548
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()
const std::vector< ReferenceCell > & get_reference_cells() const
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)
const std::vector< Point< spacedim > > & get_vertices() const
void flip_all_direction_flags()
active_cell_iterator begin_active(const unsigned int level=0) const
void execute_coarsening()
void save_refine_flags(std::vector< bool > &v) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_DEPRECATED_EARLY
Definition: config.h:167
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
Point< 2 > first
Definition: grid_out.cc:4587
unsigned int level
Definition: grid_out.cc:4590
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< active_face_iterator > active_face_iterators() const
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
IteratorRange< cell_iterator > cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
#define DeclException0(Exception0)
Definition: exceptions.h:470
static ::ExceptionBase & ExcGridReadError()
static ::ExceptionBase & ExcNonOrientableTriangulation()
static ::ExceptionBase & ExcInconsistentCoarseningFlags()
static ::ExceptionBase & ExcFacesHaveNoLevel()
static ::ExceptionBase & ExcTriangulationNotEmpty(int arg1, int arg2)
static ::ExceptionBase & ExcEmptyLevel(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:538
static ::ExceptionBase & ExcInvalidLevel(int arg1, int arg2)
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcBoundaryIdNotFound(types::boundary_id arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
static ::ExceptionBase & ExcMessage(std::string arg1)
typename IteratorSelector::hex_iterator hex_iterator
Definition: tria.h:1496
typename IteratorSelector::active_quad_iterator active_quad_iterator
Definition: tria.h:1487
typename IteratorSelector::active_hex_iterator active_hex_iterator
Definition: tria.h:1507
typename IteratorSelector::quad_iterator quad_iterator
Definition: tria.h:1472
typename IteratorSelector::line_iterator line_iterator
Definition: tria.h:1448
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition: tria.h:1361
typename IteratorSelector::active_line_iterator active_line_iterator
Definition: tria.h:1463
virtual std::vector< types::manifold_id > get_manifold_ids() const
void set_all_manifold_ids_on_boundary(const types::manifold_id number)
void reset_manifold(const types::manifold_id manifold_number)
void set_manifold(const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object)
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
void reset_all_manifolds()
void set_all_manifold_ids(const types::manifold_id number)
static const char T
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Definition: hp.h:118
Definition: types.h:32
unsigned int manifold_id
Definition: types.h:141
global_cell_index coarse_cell_id
Definition: types.h:114
unsigned int subdomain_id
Definition: types.h:43
unsigned int global_cell_index
Definition: types.h:105
unsigned int boundary_id
Definition: types.h:129
T operator()(InputIterator first, InputIterator last) const
Definition: tria.h:2064
virtual ~DistortedCellList() noexcept override
std::list< typename Triangulation< dim, spacedim >::cell_iterator > distorted_cells
Definition: tria.h:1542
boost::signals2::signal< void()> mesh_movement
Definition: tria.h:2121
boost::signals2::signal< void()> post_distributed_load
Definition: tria.h:2289
boost::signals2::signal< void()> pre_distributed_save
Definition: tria.h:2267
boost::signals2::signal< void(const Triangulation< dim, spacedim > &destination_tria)> copy
Definition: tria.h:2152
boost::signals2::signal< void()> pre_distributed_refinement
Definition: tria.h:2221
boost::signals2::signal< void()> post_distributed_refinement
Definition: tria.h:2240
boost::signals2::signal< void()> pre_refinement
Definition: tria.h:2097
boost::signals2::signal< void()> any_change
Definition: tria.h:2178
boost::signals2::signal< void()> create
Definition: tria.h:2088
boost::signals2::signal< void()> clear
Definition: tria.h:2167
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> post_refinement_on_cell
Definition: tria.h:2142
boost::signals2::signal< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell)> pre_coarsening_on_cell
Definition: tria.h:2132
boost::signals2::signal< void()> post_refinement
Definition: tria.h:2104
boost::signals2::signal< void()> pre_partition
Definition: tria.h:2112
boost::signals2::signal< void()> pre_distributed_repartition
Definition: tria.h:2252
boost::signals2::signal< void()> pre_distributed_load
Definition: tria.h:2282
boost::signals2::signal< void()> post_p4est_refinement
Definition: tria.h:2230
boost::signals2::signal< void()> post_distributed_repartition
Definition: tria.h:2259
boost::signals2::signal< unsigned int(const cell_iterator &, const CellStatus), CellWeightSum< unsigned int > > cell_weight
Definition: tria.h:2208
boost::signals2::signal< void()> post_distributed_save
Definition: tria.h:2274
void serialize(Archive &ar, const unsigned int version)
std::vector< unsigned int > n_active_lines_level
Definition: tria.h:183
std::vector< unsigned int > n_active_quads_level
Definition: tria.h:241
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:300