]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Relax assumpitions in GridTools::partition_triangulation_zorder 9103/head
authorPeter Munch <peterrmuench@gmail.com>
Wed, 27 Nov 2019 11:19:08 +0000 (12:19 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Wed, 27 Nov 2019 18:47:50 +0000 (19:47 +0100)
include/deal.II/grid/grid_tools.h
source/grid/grid_tools.cc
source/grid/grid_tools.inst.in
tests/grid/partition_triangulation_zorder_01.cc [new file with mode: 0644]
tests/grid/partition_triangulation_zorder_01.output [new file with mode: 0644]

index 6d5b276837d6eb635f26d20fbf40d48cd118db8d..d46660c28480de96d51e355d1b58ff37422b3bd0 100644 (file)
@@ -2006,15 +2006,23 @@ namespace GridTools
 
   /**
    * 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
index a328fcfe6460028a197bba820e27dc2bf69150e3..607772013f456f5dcc85de9b8763f227e637747c 100644 (file)
@@ -3037,7 +3037,8 @@ namespace GridTools
   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),
@@ -3100,40 +3101,41 @@ namespace GridTools
     // 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);
+          }
+      }
   }
 
 
index 78ff47ec8bf6a83cb39ecdf88cc64c7360ce5c19..ff06ead6babdb8422067ea21528d45c79fa96bc8 100644 (file)
@@ -297,7 +297,8 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
       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(
diff --git a/tests/grid/partition_triangulation_zorder_01.cc b/tests/grid/partition_triangulation_zorder_01.cc
new file mode 100644 (file)
index 0000000..1e1b843
--- /dev/null
@@ -0,0 +1,73 @@
+// ---------------------------------------------------------------------
+//
+// 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();
+    }
+}
diff --git a/tests/grid/partition_triangulation_zorder_01.output b/tests/grid/partition_triangulation_zorder_01.output
new file mode 100644 (file)
index 0000000..50dd0f8
--- /dev/null
@@ -0,0 +1,5 @@
+
+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 

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.