16#include <deal.II/base/mpi.templates.h>
35#ifdef DEAL_II_WITH_MPI
40 template <
int dim,
int spacedim>
54 const auto partition_settings =
63 ExcMessage(
"Settings must contain exactly one type of the active "
64 "cell partitioning scheme."));
68 ExcMessage(
"construct_multigrid_hierarchy requires "
69 "allow_artificial_cells to be set to true."));
74 template <
int dim,
int spacedim>
84 template <
int dim,
int spacedim>
92 const unsigned int max_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."));
107# ifdef DEAL_II_TRILINOS_WITH_ZOLTAN
109# elif defined DEAL_II_WITH_METIS
117# ifndef DEAL_II_TRILINOS_WITH_ZOLTAN
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."));
132# ifndef DEAL_II_WITH_METIS
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."));
185 active_halo_layer_vector =
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());
193 for (
unsigned int index = 0; cell != endc; cell++, index++)
198 if (cell->is_locally_owned() ==
false &&
199 active_halo_layer.find(cell) == active_halo_layer.end())
212 for (
unsigned int lvl = 0; lvl < this->
n_levels(); ++lvl)
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());
229 cell_iterator cell = this->
begin(lvl),
230 endc = this->
end(lvl);
231 for (
unsigned int index = 0; cell != endc; cell++, index++)
236 cell->level_subdomain_id();
245 if (cell->is_active() &&
246 cell->subdomain_id() !=
251 if (cell->has_children())
253 bool keep_cell =
false;
254 for (
unsigned int c = 0;
255 c < GeometryInfo<dim>::max_children_per_cell;
257 if (cell->child(c)->level_subdomain_id() ==
269 if (!cell->is_locally_owned_on_level() &&
270 level_halo_layer.find(cell) != level_halo_layer.
end())
274 cell->set_level_subdomain_id(
283 for (
unsigned int index = 0; cell != endc; cell++,
index++)
284 true_subdomain_ids_of_cells[
index] = cell->subdomain_id();
291 const unsigned int n_my_cells = std::count_if(
292 this->begin_active(),
295 [](
const auto &i) {
return (i.is_locally_owned()); });
297 const unsigned int total_cells =
299 Assert(total_cells == this->n_active_cells(),
300 ExcMessage(
"Not all cells are assigned to a processor."));
305 if (settings & construct_multigrid_hierarchy)
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());
313 const unsigned int total_cells =
315 Assert(total_cells == this->n_cells(),
316 ExcMessage(
"Not all cells are assigned to a processor."));
323 template <
int dim,
int spacedim>
332 template <
int dim,
int spacedim>
334 const std::vector<types::subdomain_id>
342 template <
int dim,
int spacedim>
344 const std::vector<types::subdomain_id>
346 const unsigned int level)
const
358 template <
int dim,
int spacedim>
374 using int_type = std::underlying_type_t<
376 static_assert(
sizeof(int_type) ==
sizeof(std::uint8_t),
377 "Internal type mismatch.");
379 std::vector<int_type> refinement_configurations(this->
n_active_cells() *
383 if (cell->is_locally_owned())
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);
393 refinement_configurations);
397 cell->clear_refine_flag();
398 cell->clear_coarsen_flag();
401 (refinement_configurations[cell->active_cell_index() * 2 + 0] >
405 refinement_configurations[cell->active_cell_index() * 2 +
409 "Refinement/coarsening flags of cells are not consistent in parallel!"));
411 if (refinement_configurations[cell->active_cell_index() * 2 + 0] !=
414 refinement_configurations[cell->active_cell_index() * 2 + 0]));
416 if (refinement_configurations[cell->active_cell_index() * 2 + 1] >
418 cell->set_coarsen_flag();
424 this->update_number_cache();
429 template <
int dim,
int spacedim>
442 const typename ::Triangulation<dim, spacedim>::DistortedCellList
455 template <
int dim,
int spacedim>
461 (void)construction_data;
468 template <
int dim,
int spacedim>
475 const ::parallel::DistributedTriangulationBase<dim, spacedim>
476 *
>(&other_tria) ==
nullptr),
478 "Cannot use this function on parallel::distributed::Triangulation."));
494 template <
int dim,
int spacedim>
504 template <
int dim,
int spacedim>
514 template <
int dim,
int spacedim>
516 const std::vector<unsigned int>
520 return true_subdomain_ids_of_cells;
525 template <
int dim,
int spacedim>
527 const std::vector<unsigned int>
529 const unsigned int)
const
532 return true_level_subdomain_ids_of_cells;
548 template <
int dim,
int spacedim>
560 const std::vector<types::subdomain_id> &true_subdomain_ids =
561 shared_tria->get_true_subdomain_ids_of_cells();
563 saved_subdomain_ids.resize(shared_tria->n_active_cells());
564 for (const auto &cell : shared_tria->active_cell_iterators())
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]);
575 template <
int dim,
int spacedim>
582 for (
const auto &cell :
shared_tria->active_cell_iterators())
584 const unsigned int index = cell->active_cell_index();
595#include "distributed/shared_tria.inst"
unsigned int n_active_cells() const
cell_iterator end() const
cell_iterator begin(const unsigned int level=0) const
std::vector< Point< spacedim > > vertices
unsigned int n_active_cells() const
unsigned int n_levels() const
cell_iterator end() const
unsigned int n_cells() const
MeshSmoothing smooth_grid
active_cell_iterator begin_active(const unsigned int level=0) const
std::vector< unsigned int > saved_subdomain_ids
const ObserverPointer< const ::parallel::shared::Triangulation< dim, spacedim > > shared_tria
TemporarilyRestoreSubdomainIds(const Triangulation< dim, spacedim > &tria)
types::subdomain_id n_subdomains
const MPI_Comm mpi_communicator
virtual MPI_Comm get_mpi_communicator() const override
virtual void update_number_cache()
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)
types::subdomain_id my_subdomain
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &old_tria) override
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
@ partition_custom_signal
@ construct_multigrid_hierarchy
bool with_artificial_cells() const
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
const bool allow_artificial_cells
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
constexpr bool running_in_debug_mode()
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
#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)
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