Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20: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\}}\)
Loading...
Searching...
No Matches
tria_base.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2015 - 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_distributed_tria_base_h
16#define dealii_distributed_tria_base_h
17
18
19#include <deal.II/base/config.h>
20
26
27#include <deal.II/grid/tria.h>
28
29#include <functional>
30#include <list>
31#include <set>
32#include <utility>
33#include <vector>
34
35
37
38namespace parallel
39{
77 template <int dim, int spacedim = dim>
79 class TriangulationBase : public ::Triangulation<dim, spacedim>
80 {
81 public:
86 const MPI_Comm mpi_communicator,
87 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
89 const bool check_for_distorted_cells = false);
90
94 virtual ~TriangulationBase() override;
95
99 virtual MPI_Comm
100 get_communicator() const override;
101
105 virtual bool
107
116 virtual void
117 copy_triangulation(
118 const ::Triangulation<dim, spacedim> &old_tria) override;
119
138 unsigned int
139 n_locally_owned_active_cells() const;
140
147 n_global_active_cells() const override;
148
152 virtual std::size_t
153 memory_consumption() const override;
154
155
163 virtual unsigned int
164 n_global_levels() const override;
165
173 locally_owned_subdomain() const override;
174
184 const std::set<types::subdomain_id> &
185 ghost_owners() const;
186
202 const std::set<types::subdomain_id> &
203 level_ghost_owners() const;
204
205 std::weak_ptr<const Utilities::MPI::Partitioner>
206 global_active_cell_index_partitioner() const override;
207
208 std::weak_ptr<const Utilities::MPI::Partitioner>
209 global_level_cell_index_partitioner(
210 const unsigned int level) const override;
211
227 virtual std::vector<types::boundary_id>
228 get_boundary_ids() const override;
229
245 virtual std::vector<types::manifold_id>
246 get_manifold_ids() const override;
247
300 void
301 communicate_locally_moved_vertices(
302 const std::vector<bool> &vertex_locally_moved);
303
305 n_global_coarse_cells() const override;
306
307 protected:
314
320
325
331 {
336
342
347
353 unsigned int n_global_levels;
354
359 std::set<types::subdomain_id> ghost_owners;
360
365 std::set<types::subdomain_id> level_ghost_owners;
366
370 std::shared_ptr<const Utilities::MPI::Partitioner>
372
376 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
378
379 NumberCache();
380 };
381
383
387 virtual void
388 update_number_cache();
389
390 void
391 update_reference_cells() override;
392
396 void
397 reset_global_cell_indices();
398 };
399
400
448 template <int dim, int spacedim = dim>
451 : public ::parallel::TriangulationBase<dim, spacedim>
452 {
453 public:
458 const MPI_Comm mpi_communicator,
459 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
461 const bool check_for_distorted_cells = false);
462
469 virtual void
470 clear() override;
471
473 typename ::Triangulation<dim, spacedim>::cell_iterator;
474
490 virtual bool
491 has_hanging_nodes() const override;
492
500 using Triangulation<dim, spacedim>::save;
501
502
509 using Triangulation<dim, spacedim>::load;
510
511
518 virtual void
519 load(const std::string &filename, const bool autopartition) = 0;
520 };
521
522} // namespace parallel
523
525
526#endif
virtual void load(const std::string &filename, const bool autopartition)=0
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
Definition tria_base.h:473
types::subdomain_id n_subdomains
Definition tria_base.h:324
const MPI_Comm mpi_communicator
Definition tria_base.h:313
virtual bool is_multilevel_hierarchy_constructed() const =0
types::subdomain_id my_subdomain
Definition tria_base.h:319
#define DEAL_II_DEPRECATED
Definition config.h:207
#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
unsigned int level
Definition grid_out.cc:4616
std::set< types::subdomain_id > level_ghost_owners
Definition tria_base.h:365
types::global_cell_index n_global_active_cells
Definition tria_base.h:341
std::set< types::subdomain_id > ghost_owners
Definition tria_base.h:359
std::shared_ptr< const Utilities::MPI::Partitioner > active_cell_index_partitioner
Definition tria_base.h:371
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > level_cell_index_partitioners
Definition tria_base.h:377
types::coarse_cell_id number_of_global_coarse_cells
Definition tria_base.h:346