Reference documentation for deal.II version GIT relicensing-233-g802318d791 2024-03-28 20: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
shared_tria.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_shared_tria_h
16#define dealii_distributed_shared_tria_h
17
18
19#include <deal.II/base/config.h>
20
25
27
28#include <deal.II/grid/tria.h>
29
30#include <functional>
31#include <list>
32#include <utility>
33#include <vector>
34
36
37namespace parallel
38{
39#ifdef DEAL_II_WITH_MPI
40
41
42 namespace shared
43 {
106 template <int dim, int spacedim = dim>
109 : public ::parallel::TriangulationBase<dim, spacedim>
110 {
111 public:
113 typename ::Triangulation<dim, spacedim>::active_cell_iterator;
115 typename ::Triangulation<dim, spacedim>::cell_iterator;
116
131 {
141 partition_auto = 0x0,
142
146 partition_metis = 0x1,
147
165 partition_zorder = 0x2,
166
170 partition_zoltan = 0x3,
171
221 partition_custom_signal = 0x4,
222
231 construct_multigrid_hierarchy = 0x8,
232 };
233
234
259 const MPI_Comm mpi_communicator,
260 const typename ::Triangulation<dim, spacedim>::MeshSmoothing =
262 const bool allow_artificial_cells = false,
263 const Settings settings = partition_auto);
264
268 virtual ~Triangulation() override = default;
269
273 virtual bool
274 is_multilevel_hierarchy_constructed() const override;
275
284 virtual void
285 execute_coarsening_and_refinement() override;
286
293 virtual void
294 create_triangulation(const std::vector<Point<spacedim>> &vertices,
295 const std::vector<CellData<dim>> &cells,
296 const SubCellData &subcelldata) override;
297
303 virtual void
304 create_triangulation(
306 &construction_data) override;
307
318 virtual void
319 copy_triangulation(
320 const ::Triangulation<dim, spacedim> &other_tria) override;
321
329 using parallel::TriangulationBase<dim, spacedim>::save;
330
336 using parallel::TriangulationBase<dim, spacedim>::load;
337
346 template <class Archive>
347 void
348 load(Archive &ar, const unsigned int version);
349
358 const std::vector<types::subdomain_id> &
359 get_true_subdomain_ids_of_cells() const;
360
369 const std::vector<types::subdomain_id> &
370 get_true_level_subdomain_ids_of_cells(const unsigned int level) const;
371
376 bool
377 with_artificial_cells() const;
378
379 private:
384
389
394 void
395 partition();
396
409 std::vector<types::subdomain_id> true_subdomain_ids_of_cells;
410
419 std::vector<std::vector<types::subdomain_id>>
421 };
422
423
424
425 template <int dim, int spacedim>
427 template <class Archive>
429 const unsigned int version)
430 {
432 partition();
433 this->update_number_cache();
434 }
435 } // namespace shared
436#else
437
438 namespace shared
439 {
453 template <int dim, int spacedim = dim>
455 class Triangulation
456 : public ::parallel::TriangulationBase<dim, spacedim>
457 {
458 public:
463 Triangulation() = delete;
464
468 virtual bool
470
474 const std::vector<types::subdomain_id> &
476
480 const std::vector<types::subdomain_id> &
481 get_true_level_subdomain_ids_of_cells(const unsigned int level) const;
482
486 bool
487 with_artificial_cells() const;
488
489 private:
493 std::vector<types::subdomain_id> true_subdomain_ids_of_cells;
494
498 std::vector<types::subdomain_id> true_level_subdomain_ids_of_cells;
499 };
500 } // namespace shared
501
502
503#endif
504} // namespace parallel
505
506
507namespace internal
508{
509 namespace parallel
510 {
511 namespace shared
512 {
532 template <int dim, int spacedim = dim>
534 {
535 public:
545 const Triangulation<dim, spacedim> &tria);
546
554
555 private:
559 const SmartPointer<
560 const ::parallel::shared::Triangulation<dim, spacedim>>
562
568 std::vector<unsigned int> saved_subdomain_ids;
569 };
570 } // namespace shared
571 } // namespace parallel
572} // namespace internal
573
575
576#endif
Definition point.h:111
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:502
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
Point< 3 > vertices[4]
unsigned int level
Definition grid_out.cc:4616