deal.II version 9.7.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.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2015 - 2025 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#include <deal.II/base/mpi.h>
16#include <deal.II/base/mpi.templates.h>
18
21
24#include <deal.II/grid/tria.h>
27
29
30#include <type_traits>
31
32
34
35#ifdef DEAL_II_WITH_MPI
36namespace parallel
37{
38 namespace shared
39 {
40 template <int dim, int spacedim>
41 DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
44 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
46 const bool allow_artificial_cells,
47 const Settings settings)
50 false)
53 {
54 const auto partition_settings =
58 Assert(partition_settings == partition_auto ||
59 partition_settings == partition_metis ||
60 partition_settings == partition_zoltan ||
61 partition_settings == partition_zorder ||
62 partition_settings == partition_custom_signal,
63 ExcMessage("Settings must contain exactly one type of the active "
64 "cell partitioning scheme."));
65
68 ExcMessage("construct_multigrid_hierarchy requires "
69 "allow_artificial_cells to be set to true."));
70 }
71
72
73
74 template <int dim, int spacedim>
81
82
83
84 template <int dim, int spacedim>
86 void Triangulation<dim, spacedim>::partition()
87 {
88 if constexpr (running_in_debug_mode())
89 {
90 // Check that all meshes are the same (or at least have the same
91 // total number of active cells):
92 const unsigned int max_active_cells =
94 this->get_mpi_communicator());
95 Assert(
96 max_active_cells == this->n_active_cells(),
98 "A parallel::shared::Triangulation needs to be refined in the same "
99 "way on all processors, but the participating processors don't "
100 "agree on the number of active cells."));
101 }
102
103 auto partition_settings = (partition_zoltan | partition_metis |
105 settings;
106 if (partition_settings == partition_auto)
107# ifdef DEAL_II_TRILINOS_WITH_ZOLTAN
108 partition_settings = partition_zoltan;
109# elif defined DEAL_II_WITH_METIS
110 partition_settings = partition_metis;
111# else
112 partition_settings = partition_zorder;
113# endif
114
115 if (partition_settings == partition_zoltan)
116 {
117# ifndef DEAL_II_TRILINOS_WITH_ZOLTAN
118 AssertThrow(false,
120 "Choosing 'partition_zoltan' requires the library "
121 "to be compiled with support for Zoltan! "
122 "Instead, you might use 'partition_auto' to select "
123 "a partitioning algorithm that is supported "
124 "by your current configuration."));
125# else
128# endif
129 }
130 else if (partition_settings == partition_metis)
131 {
132# ifndef DEAL_II_WITH_METIS
133 AssertThrow(false,
135 "Choosing 'partition_metis' requires the library "
136 "to be compiled with support for METIS! "
137 "Instead, you might use 'partition_auto' to select "
138 "a partitioning algorithm that is supported "
139 "by your current configuration."));
140# else
142 *this,
144# endif
145 }
146 else if (partition_settings == partition_zorder)
147 {
149 }
150 else if (partition_settings == partition_custom_signal)
151 {
152 // User partitions mesh manually
153 }
154 else
155 {
157 }
158
159 // do not partition multigrid levels if user is
160 // defining a custom partition
164
166
167 // loop over all cells and mark artificial:
170 cell = this->begin_active(),
171 endc = this->end();
172
174 {
175 // get active halo layer of (ghost) cells
176 // parallel::shared::Triangulation<dim>::
177 std::function<bool(
181
182 const std::vector<typename parallel::shared::Triangulation<
183 dim,
184 spacedim>::active_cell_iterator>
185 active_halo_layer_vector =
187 predicate);
188 std::set<typename parallel::shared::Triangulation<dim, spacedim>::
189 active_cell_iterator>
190 active_halo_layer(active_halo_layer_vector.begin(),
191 active_halo_layer_vector.end());
192
193 for (unsigned int index = 0; cell != endc; cell++, index++)
194 {
195 // store original/true subdomain ids:
196 true_subdomain_ids_of_cells[index] = cell->subdomain_id();
197
198 if (cell->is_locally_owned() == false &&
199 active_halo_layer.find(cell) == active_halo_layer.end())
200 cell->set_subdomain_id(numbers::artificial_subdomain_id);
201 }
202
203 // loop over all cells in multigrid hierarchy and mark artificial:
205 {
207
208 std::function<bool(
210 cell_iterator &)>
212 for (unsigned int lvl = 0; lvl < this->n_levels(); ++lvl)
213 {
215 this->n_cells(lvl));
216
217 const std::vector<typename parallel::shared::Triangulation<
218 dim,
219 spacedim>::cell_iterator>
220 level_halo_layer_vector =
222 *this, predicate, lvl);
223 std::set<typename parallel::shared::
224 Triangulation<dim, spacedim>::cell_iterator>
225 level_halo_layer(level_halo_layer_vector.begin(),
226 level_halo_layer_vector.end());
227
229 cell_iterator cell = this->begin(lvl),
230 endc = this->end(lvl);
231 for (unsigned int index = 0; cell != endc; cell++, index++)
232 {
233 // Store true level subdomain IDs before setting
234 // artificial
236 cell->level_subdomain_id();
237
238 // for active cells, we must have knowledge of level
239 // subdomain ids of all neighbors to our subdomain, not
240 // just neighbors on the same level. if the cells
241 // subdomain id was not set to artitficial above, we will
242 // also keep its level subdomain id since it is either
243 // owned by this processor or in the ghost layer of the
244 // active mesh.
245 if (cell->is_active() &&
246 cell->subdomain_id() !=
248 continue;
249
250 // we must have knowledge of our parent in the hierarchy
251 if (cell->has_children())
252 {
253 bool keep_cell = false;
254 for (unsigned int c = 0;
255 c < GeometryInfo<dim>::max_children_per_cell;
256 ++c)
257 if (cell->child(c)->level_subdomain_id() ==
259 {
260 keep_cell = true;
261 break;
262 }
263 if (keep_cell)
264 continue;
265 }
266
267 // we must have knowledge of our neighbors on the same
268 // level
269 if (!cell->is_locally_owned_on_level() &&
270 level_halo_layer.find(cell) != level_halo_layer.end())
271 continue;
272
273 // mark all other cells to artificial
274 cell->set_level_subdomain_id(
276 }
277 }
278 }
279 }
280 else
281 {
282 // just store true subdomain ids
283 for (unsigned int index = 0; cell != endc; cell++, index++)
284 true_subdomain_ids_of_cells[index] = cell->subdomain_id();
286
287 if constexpr (running_in_debug_mode())
288 {
289 {
290 // Assert that each cell is owned by a processor
291 const unsigned int n_my_cells = std::count_if(
292 this->begin_active(),
294 this->end()),
295 [](const auto &i) { return (i.is_locally_owned()); });
296
297 const unsigned int total_cells =
298 Utilities::MPI::sum(n_my_cells, this->get_mpi_communicator());
299 Assert(total_cells == this->n_active_cells(),
300 ExcMessage("Not all cells are assigned to a processor."));
301 }
302
303 // If running with multigrid, assert that each level
304 // cell is owned by a processor
305 if (settings & construct_multigrid_hierarchy)
306 {
307 const unsigned int n_my_cells =
308 std::count_if(this->begin(), this->end(), [](const auto &i) {
309 return (i.is_locally_owned_on_level());
310 });
311
312
313 const unsigned int total_cells =
314 Utilities::MPI::sum(n_my_cells, this->get_mpi_communicator());
315 Assert(total_cells == this->n_cells(),
316 ExcMessage("Not all cells are assigned to a processor."));
317 }
318 }
320
321
322
323 template <int dim, int spacedim>
326 {
328 }
329
330
331
332 template <int dim, int spacedim>
334 const std::vector<types::subdomain_id>
339
340
341
342 template <int dim, int spacedim>
344 const std::vector<types::subdomain_id>
346 const unsigned int level) const
347 {
351 this->n_cells(level),
354 }
355
356
357
358 template <int dim, int spacedim>
361 {
362 // make sure that all refinement/coarsening flags are the same on all
363 // processes
364 {
365 // Obtain the type used to store the different possibilities
366 // a cell can be refined. This is a bit awkward because
367 // what `cell->refine_flag_set()` returns is a struct
368 // type, RefinementCase, which internally stores a
369 // std::uint8_t, which actually holds integers of
370 // enum type RefinementPossibilities<dim>::Possibilities.
371 // In the following, use the actual name of the enum, but
372 // make sure that it is in fact a `std::uint8_t` or
373 // equally sized type.
374 using int_type = std::underlying_type_t<
376 static_assert(sizeof(int_type) == sizeof(std::uint8_t),
377 "Internal type mismatch.");
378
379 std::vector<int_type> refinement_configurations(this->n_active_cells() *
380 2,
381 int_type(0));
382 for (const auto &cell : this->active_cell_iterators())
383 if (cell->is_locally_owned())
384 {
385 refinement_configurations[cell->active_cell_index() * 2 + 0] =
386 static_cast<int_type>(cell->refine_flag_set());
387 refinement_configurations[cell->active_cell_index() * 2 + 1] =
388 static_cast<int_type>(cell->coarsen_flag_set() ? 1 : 0);
389 }
390
391 Utilities::MPI::max(refinement_configurations,
393 refinement_configurations);
394
395 for (const auto &cell : this->active_cell_iterators())
396 {
397 cell->clear_refine_flag();
398 cell->clear_coarsen_flag();
399
400 Assert(
401 (refinement_configurations[cell->active_cell_index() * 2 + 0] >
402 0 ?
403 1 :
404 0) +
405 refinement_configurations[cell->active_cell_index() * 2 +
406 1] <=
407 1,
409 "Refinement/coarsening flags of cells are not consistent in parallel!"));
411 if (refinement_configurations[cell->active_cell_index() * 2 + 0] !=
412 0)
413 cell->set_refine_flag(RefinementCase<dim>(
414 refinement_configurations[cell->active_cell_index() * 2 + 0]));
415
416 if (refinement_configurations[cell->active_cell_index() * 2 + 1] >
417 0)
418 cell->set_coarsen_flag();
419 }
420 }
421
423 partition();
424 this->update_number_cache();
425 }
426
427
428
429 template <int dim, int spacedim>
432 const std::vector<Point<spacedim>> &vertices,
433 const std::vector<CellData<dim>> &cells,
434 const SubCellData &subcelldata)
435 {
436 try
437 {
439 vertices, cells, subcelldata);
440 }
441 catch (
442 const typename ::Triangulation<dim, spacedim>::DistortedCellList
443 &)
444 {
445 // the underlying triangulation should not be checking for distorted
446 // cells
448 }
449 partition();
450 this->update_number_cache();
451 }
452
453
454
455 template <int dim, int spacedim>
458 const TriangulationDescription::Description<dim, spacedim>
459 &construction_data)
460 {
461 (void)construction_data;
462
464 }
465
466
467
468 template <int dim, int spacedim>
470 void Triangulation<dim, spacedim>::copy_triangulation(
471 const ::Triangulation<dim, spacedim> &other_tria)
472 {
473 Assert(
474 (dynamic_cast<
475 const ::parallel::DistributedTriangulationBase<dim, spacedim>
476 *>(&other_tria) == nullptr),
478 "Cannot use this function on parallel::distributed::Triangulation."));
479
481 other_tria);
482 partition();
483 this->update_number_cache();
484 }
485 } // namespace shared
486} // namespace parallel
487
488#else
489
490namespace parallel
491{
492 namespace shared
493 {
494 template <int dim, int spacedim>
497 {
499 return true;
500 }
501
502
503
504 template <int dim, int spacedim>
505 DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
507 const
508 {
509 return false;
510 }
511
512
513
514 template <int dim, int spacedim>
515 DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
516 const std::vector<unsigned int>
518 {
520 return true_subdomain_ids_of_cells;
521 }
522
523
524
525 template <int dim, int spacedim>
526 DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
527 const std::vector<unsigned int>
529 const unsigned int) const
530 {
532 return true_level_subdomain_ids_of_cells;
533 }
534 } // namespace shared
535} // namespace parallel
536
537
538#endif
539
540
541
542namespace internal
543{
544 namespace parallel
545 {
546 namespace shared
547 {
548 template <int dim, int spacedim>
551 : shared_tria(
552 dynamic_cast<
553 const ::parallel::shared::Triangulation<dim, spacedim> *>(
554 &tria))
555 {
556 if (shared_tria && shared_tria->with_artificial_cells())
557 {
558 // Save the current set of subdomain IDs, and set subdomain IDs
559 // to the "true" owner of each cell.
560 const std::vector<types::subdomain_id> &true_subdomain_ids =
561 shared_tria->get_true_subdomain_ids_of_cells();
562
563 saved_subdomain_ids.resize(shared_tria->n_active_cells());
564 for (const auto &cell : shared_tria->active_cell_iterators())
565 {
566 const unsigned int index = cell->active_cell_index();
567 saved_subdomain_ids[index] = cell->subdomain_id();
568 cell->set_subdomain_id(true_subdomain_ids[index]);
569 }
570 }
571 }
572
573
574
575 template <int dim, int spacedim>
578 {
579 if (shared_tria && shared_tria->with_artificial_cells())
580 {
581 // Undo the subdomain modification.
582 for (const auto &cell : shared_tria->active_cell_iterators())
583 {
584 const unsigned int index = cell->active_cell_index();
585 cell->set_subdomain_id(saved_subdomain_ids[index]);
586 }
587 }
588 }
589 } // namespace shared
590 } // namespace parallel
591} // namespace internal
592
593
594/*-------------- Explicit Instantiations -------------------------------*/
595#include "distributed/shared_tria.inst"
596
Definition point.h:113
cell_iterator end() const
cell_iterator begin(const unsigned int level=0) const
std::vector< Point< spacedim > > vertices
Definition tria.h:4498
active_cell_iterator begin_active(const unsigned int level=0) const
const ObserverPointer< const ::parallel::shared::Triangulation< dim, spacedim > > shared_tria
TemporarilyRestoreSubdomainIds(const Triangulation< dim, spacedim > &tria)
virtual MPI_Comm get_mpi_communicator() const override
Definition tria_base.cc:160
TriangulationBase(const MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const bool check_for_distorted_cells=false)
Definition tria_base.cc:47
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &old_tria) override
Definition tria_base.cc:67
virtual void execute_coarsening_and_refinement() override
std::vector< types::subdomain_id > true_subdomain_ids_of_cells
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &other_tria) override
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
const std::vector< types::subdomain_id > & get_true_level_subdomain_ids_of_cells(const unsigned int level) const
Triangulation(const MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing=(::Triangulation< dim, spacedim >::none), const bool allow_artificial_cells=false, const Settings settings=partition_auto)
typename ::Triangulation< dim, spacedim >::active_cell_iterator active_cell_iterator
virtual bool is_multilevel_hierarchy_constructed() const override
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata) override
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
constexpr bool running_in_debug_mode()
Definition config.h:78
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:248
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
IteratorRange< active_cell_iterator > active_cell_iterators() const
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void partition_triangulation_zorder(const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const bool group_siblings=true)
std::vector< typename MeshType::active_cell_iterator > compute_active_cell_halo_layer(const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate)
void partition_multigrid_levels(Triangulation< dim, spacedim > &triangulation)
std::vector< typename MeshType::cell_iterator > compute_cell_halo_layer_on_level(const MeshType &mesh, const std::function< bool(const typename MeshType::cell_iterator &)> &predicate, const unsigned int level)
void partition_triangulation(const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner=SparsityTools::Partitioner::metis)
T sum(const T &t, const MPI_Comm mpi_communicator)
T max(const T &t, const MPI_Comm mpi_communicator)
constexpr types::subdomain_id artificial_subdomain_id
Definition types.h:402
STL namespace.