Reference documentation for deal.II version GIT f6a5d312c9 2023-10-04 08:50:02+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\}}\)
tria_base.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2023 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_distributed_tria_base_h
17 #define dealii_distributed_tria_base_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/mpi_stub.h>
27 
28 #include <deal.II/grid/tria.h>
29 
30 #include <functional>
31 #include <list>
32 #include <set>
33 #include <utility>
34 #include <vector>
35 
36 
38 
39 namespace parallel
40 {
78  template <int dim, int spacedim = dim>
79  DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
80  class TriangulationBase : public ::Triangulation<dim, spacedim>
81  {
82  public:
87  const MPI_Comm mpi_communicator,
88  const typename ::Triangulation<dim, spacedim>::MeshSmoothing
89  smooth_grid = (::Triangulation<dim, spacedim>::none),
90  const bool check_for_distorted_cells = false);
91 
95  virtual ~TriangulationBase() override;
96 
100  virtual MPI_Comm
101  get_communicator() const override;
102 
106  virtual bool
108 
117  virtual void
118  copy_triangulation(
119  const ::Triangulation<dim, spacedim> &old_tria) override;
120 
139  unsigned int
140  n_locally_owned_active_cells() const;
141 
148  n_global_active_cells() const override;
149 
153  virtual std::size_t
154  memory_consumption() const override;
155 
156 
164  virtual unsigned int
165  n_global_levels() const override;
166 
174  locally_owned_subdomain() const override;
175 
185  const std::set<types::subdomain_id> &
186  ghost_owners() const;
187 
203  const std::set<types::subdomain_id> &
204  level_ghost_owners() const;
205 
206  std::weak_ptr<const Utilities::MPI::Partitioner>
207  global_active_cell_index_partitioner() const override;
208 
209  std::weak_ptr<const Utilities::MPI::Partitioner>
210  global_level_cell_index_partitioner(
211  const unsigned int level) const override;
212 
219  virtual std::vector<types::boundary_id>
220  get_boundary_ids() const override;
221 
228  virtual std::vector<types::manifold_id>
229  get_manifold_ids() const override;
230 
283  void
284  communicate_locally_moved_vertices(
285  const std::vector<bool> &vertex_locally_moved);
286 
287  virtual types::coarse_cell_id
288  n_global_coarse_cells() const override;
289 
290  protected:
296  const MPI_Comm mpi_communicator;
297 
303 
308 
313  struct NumberCache
314  {
319 
325 
330 
336  unsigned int n_global_levels;
337 
342  std::set<types::subdomain_id> ghost_owners;
343 
348  std::set<types::subdomain_id> level_ghost_owners;
349 
353  std::shared_ptr<const Utilities::MPI::Partitioner>
355 
359  std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
361 
362  NumberCache();
363  };
364 
366 
370  virtual void
371  update_number_cache();
372 
376  void
377  update_reference_cells() override;
378 
382  void
383  reset_global_cell_indices();
384  };
385 
386 
434  template <int dim, int spacedim = dim>
435  DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
437  : public ::parallel::TriangulationBase<dim, spacedim>
438  {
439  public:
444  const MPI_Comm mpi_communicator,
445  const typename ::Triangulation<dim, spacedim>::MeshSmoothing
446  smooth_grid = (::Triangulation<dim, spacedim>::none),
447  const bool check_for_distorted_cells = false);
448 
455  virtual void
456  clear() override;
457 
459  typename ::Triangulation<dim, spacedim>::cell_iterator;
460 
476  virtual bool
477  has_hanging_nodes() const override;
478 
487 
488 
496 
497 
504  virtual void
505  load(const std::string &filename, const bool autopartition) = 0;
506  };
507 
508 } // namespace parallel
509 
511 
512 #endif
virtual void load(const std::string &filename, const bool autopartition)=0
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
Definition: tria_base.h:459
types::subdomain_id n_subdomains
Definition: tria_base.h:307
const MPI_Comm mpi_communicator
Definition: tria_base.h:296
virtual bool is_multilevel_hierarchy_constructed() const =0
types::subdomain_id my_subdomain
Definition: tria_base.h:302
#define DEAL_II_DEPRECATED
Definition: config.h:178
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:477
#define DEAL_II_CXX20_REQUIRES(condition)
Definition: config.h:166
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:478
unsigned int level
Definition: grid_out.cc:4617
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
global_cell_index coarse_cell_id
Definition: types.h:127
unsigned int subdomain_id
Definition: types.h:44
unsigned int global_cell_index
Definition: types.h:118
std::set< types::subdomain_id > level_ghost_owners
Definition: tria_base.h:348
types::global_cell_index n_global_active_cells
Definition: tria_base.h:324
std::set< types::subdomain_id > ghost_owners
Definition: tria_base.h:342
std::shared_ptr< const Utilities::MPI::Partitioner > active_cell_index_partitioner
Definition: tria_base.h:354
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > level_cell_index_partitioners
Definition: tria_base.h:360
types::coarse_cell_id number_of_global_coarse_cells
Definition: tria_base.h:329