Reference documentation for deal.II version 9.5.0
\(\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
shared_tria.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_shared_tria_h
17#define dealii_distributed_shared_tria_h
18
19
20#include <deal.II/base/config.h>
21
26
28
29#include <deal.II/grid/tria.h>
30
31#include <functional>
32#include <list>
33#include <utility>
34#include <vector>
35
37
38namespace parallel
39{
40#ifdef DEAL_II_WITH_MPI
41
42
43 namespace shared
44 {
107 template <int dim, int spacedim = dim>
110 : public ::parallel::TriangulationBase<dim, spacedim>
111 {
112 public:
114 typename ::Triangulation<dim, spacedim>::active_cell_iterator;
116 typename ::Triangulation<dim, spacedim>::cell_iterator;
117
132 {
142 partition_auto = 0x0,
143
147 partition_metis = 0x1,
148
166 partition_zorder = 0x2,
167
171 partition_zoltan = 0x3,
172
222 partition_custom_signal = 0x4,
223
232 construct_multigrid_hierarchy = 0x8,
233 };
234
235
260 const MPI_Comm mpi_communicator,
261 const typename ::Triangulation<dim, spacedim>::MeshSmoothing =
263 const bool allow_artificial_cells = false,
264 const Settings settings = partition_auto);
265
269 virtual ~Triangulation() override = default;
270
274 virtual bool
275 is_multilevel_hierarchy_constructed() const override;
276
285 virtual void
286 execute_coarsening_and_refinement() override;
287
294 virtual void
295 create_triangulation(const std::vector<Point<spacedim>> &vertices,
296 const std::vector<CellData<dim>> & cells,
297 const SubCellData &subcelldata) override;
298
304 virtual void
305 create_triangulation(
307 &construction_data) override;
308
319 virtual void
320 copy_triangulation(
321 const ::Triangulation<dim, spacedim> &other_tria) override;
322
331 template <class Archive>
332 void
333 load(Archive &ar, const unsigned int version);
334
343 const std::vector<types::subdomain_id> &
344 get_true_subdomain_ids_of_cells() const;
345
354 const std::vector<types::subdomain_id> &
355 get_true_level_subdomain_ids_of_cells(const unsigned int level) const;
356
361 bool
362 with_artificial_cells() const;
363
364 private:
369
374
379 void
380 partition();
381
394 std::vector<types::subdomain_id> true_subdomain_ids_of_cells;
395
404 std::vector<std::vector<types::subdomain_id>>
406 };
407
408
409
410 template <int dim, int spacedim>
412 template <class Archive>
414 const unsigned int version)
415 {
417 partition();
418 this->update_number_cache();
419 }
420 } // namespace shared
421#else
422
423 namespace shared
424 {
438 template <int dim, int spacedim = dim>
440 class Triangulation
441 : public ::parallel::TriangulationBase<dim, spacedim>
442 {
443 public:
448 Triangulation() = delete;
449
453 virtual bool
455
459 const std::vector<types::subdomain_id> &
461
465 const std::vector<types::subdomain_id> &
466 get_true_level_subdomain_ids_of_cells(const unsigned int level) const;
467
471 bool
472 with_artificial_cells() const;
473
474 private:
478 std::vector<types::subdomain_id> true_subdomain_ids_of_cells;
479
483 std::vector<types::subdomain_id> true_level_subdomain_ids_of_cells;
484 };
485 } // namespace shared
486
487
488#endif
489} // namespace parallel
490
491
492namespace internal
493{
494 namespace parallel
495 {
496 namespace shared
497 {
517 template <int dim, int spacedim = dim>
519 {
520 public:
531
539
540 private:
544 const SmartPointer<
545 const ::parallel::shared::Triangulation<dim, spacedim>>
547
553 std::vector<unsigned int> saved_subdomain_ids;
554 };
555 } // namespace shared
556 } // namespace parallel
557} // namespace internal
558
560
561#endif
Definition point.h:112
const SmartPointer< const ::parallel::shared::Triangulation< dim, spacedim > > shared_tria
std::vector< types::subdomain_id > true_subdomain_ids_of_cells
typename ::Triangulation< dim, spacedim >::active_cell_iterator active_cell_iterator
const std::vector< types::subdomain_id > & get_true_subdomain_ids_of_cells() const
std::vector< std::vector< types::subdomain_id > > true_level_subdomain_ids_of_cells
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
const std::vector< types::subdomain_id > & get_true_level_subdomain_ids_of_cells(const unsigned int level) const
virtual ~Triangulation() override=default
virtual bool is_multilevel_hierarchy_constructed() const override
void load(Archive &ar, const unsigned int version)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:160
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
Point< 3 > vertices[4]
unsigned int level
Definition grid_out.cc:4618
const TriangulationDescription::Settings settings
const ::Triangulation< dim, spacedim > & tria