/**
* Generates a partitioning of the active cells making up the entire domain
- * using the same partitioning scheme as in the p4est library. After calling
- * this function, the subdomain ids of all active cells will have values
- * between zero and @p n_partitions-1. You can access the subdomain id of a
- * cell by using <tt>cell-@>subdomain_id()</tt>.
+ * using the same partitioning scheme as in the p4est library if the flag
+ * @p group_siblings is set to true (default behavior of this function).
+ * After calling this function, the subdomain ids of all active cells will
+ * have values between zero and @p n_partitions-1. You can access the
+ * subdomain id of a cell by using <tt>cell-@>subdomain_id()</tt>.
+ *
+ * @note If the flag @p group_siblings is set to false, children of a
+ * cell might be placed on different processors even though they are all
+ * active, which is an assumption made by p4est. By relaxing this, we
+ * can can create partitions owning a single cell (also for refined
+ * meshes).
*/
template <int dim, int spacedim>
void
partition_triangulation_zorder(const unsigned int n_partitions,
- Triangulation<dim, spacedim> &triangulation);
+ Triangulation<dim, spacedim> &triangulation,
+ const bool group_siblings = true);
/**
* Partitions the cells of a multigrid hierarchy by assigning level subdomain
template <int dim, int spacedim>
void
partition_triangulation_zorder(const unsigned int n_partitions,
- Triangulation<dim, spacedim> &triangulation)
+ Triangulation<dim, spacedim> &triangulation,
+ const bool group_siblings)
{
Assert((dynamic_cast<parallel::distributed::Triangulation<dim, spacedim> *>(
&triangulation) == nullptr),
// the processor with the largest number of children
// (ties are broken by picking the lower rank).
// Duplicate this logic here.
- {
- typename Triangulation<dim, spacedim>::cell_iterator
- cell = triangulation.begin(),
- endc = triangulation.end();
- for (; cell != endc; ++cell)
- {
- if (cell->active())
- continue;
- bool all_children_active = true;
- std::map<unsigned int, unsigned int> map_cpu_n_cells;
- for (unsigned int n = 0; n < cell->n_children(); ++n)
- if (!cell->child(n)->active())
- {
- all_children_active = false;
- break;
- }
- else
- ++map_cpu_n_cells[cell->child(n)->subdomain_id()];
+ if (group_siblings)
+ {
+ typename Triangulation<dim, spacedim>::cell_iterator
+ cell = triangulation.begin(),
+ endc = triangulation.end();
+ for (; cell != endc; ++cell)
+ {
+ if (cell->active())
+ continue;
+ bool all_children_active = true;
+ std::map<unsigned int, unsigned int> map_cpu_n_cells;
+ for (unsigned int n = 0; n < cell->n_children(); ++n)
+ if (!cell->child(n)->active())
+ {
+ all_children_active = false;
+ break;
+ }
+ else
+ ++map_cpu_n_cells[cell->child(n)->subdomain_id()];
- if (!all_children_active)
- continue;
+ if (!all_children_active)
+ continue;
- unsigned int new_owner = cell->child(0)->subdomain_id();
- for (std::map<unsigned int, unsigned int>::iterator it =
- map_cpu_n_cells.begin();
- it != map_cpu_n_cells.end();
- ++it)
- if (it->second > map_cpu_n_cells[new_owner])
- new_owner = it->first;
+ unsigned int new_owner = cell->child(0)->subdomain_id();
+ for (std::map<unsigned int, unsigned int>::iterator it =
+ map_cpu_n_cells.begin();
+ it != map_cpu_n_cells.end();
+ ++it)
+ if (it->second > map_cpu_n_cells[new_owner])
+ new_owner = it->first;
- for (unsigned int n = 0; n < cell->n_children(); ++n)
- cell->child(n)->set_subdomain_id(new_owner);
- }
- }
+ for (unsigned int n = 0; n < cell->n_children(); ++n)
+ cell->child(n)->set_subdomain_id(new_owner);
+ }
+ }
}
template void
partition_triangulation_zorder(
const unsigned int,
- Triangulation<deal_II_dimension, deal_II_space_dimension> &);
+ Triangulation<deal_II_dimension, deal_II_space_dimension> &,
+ const bool);
template void
partition_multigrid_levels(
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// Test GridTools::partition_triangulation_zorder
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test(const int n_refinements, const int n_partitions, const bool blocked)
+{
+ // create serial triangulation
+ Triangulation<dim> tria;
+ GridGenerator::hyper_cube(tria);
+ tria.refine_global(n_refinements);
+
+ GridTools::partition_triangulation_zorder(n_partitions, tria, blocked);
+
+ for (const auto &cell : tria.active_cell_iterators())
+ deallog << cell->subdomain_id() << " ";
+ deallog << std::endl;
+}
+
+int
+main(int argc, char *argv[])
+{
+ initlog();
+
+ if (true)
+ {
+ deallog.push("2d");
+ test<2>(1, 4, false);
+ deallog.pop();
+ }
+ if (true)
+ {
+ deallog.push("2d-pdt");
+ test<2>(1, 4, true);
+ deallog.pop();
+ }
+ if (true)
+ {
+ deallog.push("3d");
+ test<3>(1, 8, false);
+ deallog.pop();
+ }
+ if (true)
+ {
+ deallog.push("3d-pdt");
+ test<3>(1, 8, true);
+ deallog.pop();
+ }
+}
--- /dev/null
+
+DEAL:2d::0 1 2 3
+DEAL:2d-pdt::0 0 0 0
+DEAL:3d::0 1 2 3 4 5 6 7
+DEAL:3d-pdt::0 0 0 0 0 0 0 0