]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Extend RepartitioningPolicyTools::DefaultPolicy 12810/head
authorPeter Munch <peterrmuench@gmail.com>
Sat, 9 Oct 2021 17:40:44 +0000 (19:40 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Sat, 9 Oct 2021 19:06:13 +0000 (21:06 +0200)
include/deal.II/distributed/repartitioning_policy_tools.h
source/distributed/repartitioning_policy_tools.cc
tests/fullydistributed_grids/repartitioning_07.cc [new file with mode: 0644]
tests/fullydistributed_grids/repartitioning_07.with_p4est=true.mpirun=4.output [new file with mode: 0644]

index 3bdcf779455830dddea62dfb4a09d5bfb5af7932..ff283c56fdbf978f5db5dc18350dbcec6de78970 100644 (file)
@@ -67,9 +67,20 @@ namespace RepartitioningPolicyTools
   class DefaultPolicy : public Base<dim, spacedim>
   {
   public:
+    /**
+     * Constructor.
+     *
+     * @param tighten allows to renumber of subdomains so that empty ranks are
+     *   positioned at the end.
+     */
+    DefaultPolicy(const bool tighten = false);
+
     virtual LinearAlgebra::distributed::Vector<double>
     partition(
       const Triangulation<dim, spacedim> &tria_coarse_in) const override;
+
+  private:
+    const bool tighten;
   };
 
   /**
index 92d39e215ca3b3a0e9b1010c1cc579879d64b858..9b203f93c2abeb9d067b12dc76aac4051d64573f 100644 (file)
@@ -46,12 +46,61 @@ namespace RepartitioningPolicyTools
     }
   } // namespace
 
+
+  template <int dim, int spacedim>
+  DefaultPolicy<dim, spacedim>::DefaultPolicy(const bool tighten)
+    : tighten(tighten)
+  {}
+
   template <int dim, int spacedim>
   LinearAlgebra::distributed::Vector<double>
   DefaultPolicy<dim, spacedim>::partition(
-    const Triangulation<dim, spacedim> &) const
+    const Triangulation<dim, spacedim> &tria_in) const
   {
-    return {}; // nothing to do
+    if (tighten == false)
+      return {}; // nothing to do
+
+    const auto tria =
+      dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
+        &tria_in);
+
+    if (tria == nullptr)
+      return {}; // nothing to do, since serial triangulation
+
+#ifndef DEAL_II_WITH_MPI
+    Assert(false, ExcNeedsMPI());
+    return {};
+#else
+
+    const auto comm = tria->get_communicator();
+
+    const unsigned int process_has_active_locally_owned_cells =
+      tria->n_locally_owned_active_cells() > 0;
+    const unsigned int n_processes_with_active_locally_owned_cells =
+      Utilities::MPI::sum(process_has_active_locally_owned_cells, comm);
+
+    if (n_processes_with_active_locally_owned_cells ==
+        Utilities::MPI::n_mpi_processes(comm))
+      return {}; // nothing to do, since all processes have cells
+
+    unsigned int offset = 0;
+
+    const int ierr = MPI_Exscan(&process_has_active_locally_owned_cells,
+                                &offset,
+                                1,
+                                Utilities::MPI::internal::mpi_type_id(
+                                  &process_has_active_locally_owned_cells),
+                                MPI_SUM,
+                                comm);
+    AssertThrowMPI(ierr);
+
+    LinearAlgebra::distributed::Vector<double> partition(
+      tria->global_active_cell_index_partitioner().lock());
+
+    partition = offset;
+
+    return partition;
+#endif
   }
 
 
diff --git a/tests/fullydistributed_grids/repartitioning_07.cc b/tests/fullydistributed_grids/repartitioning_07.cc
new file mode 100644 (file)
index 0000000..091107a
--- /dev/null
@@ -0,0 +1,96 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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 RepartitioningPolicyTools::DefaultPolicy's tightening functionality.
+
+#include <deal.II/base/mpi_consensus_algorithms.h>
+
+#include <deal.II/distributed/fully_distributed_tria.h>
+#include <deal.II/distributed/repartitioning_policy_tools.h>
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#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_description.h>
+
+#include <deal.II/lac/la_parallel_vector.h>
+
+#include <deal.II/numerics/data_out.h>
+
+#include "./tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test(const MPI_Comm comm)
+{
+  parallel::shared::Triangulation<dim> tria(
+    comm,
+    Triangulation<dim>::none,
+    true,
+    parallel::shared::Triangulation<dim>::partition_custom_signal);
+
+  tria.signals.create.connect([&tria]() {
+    for (const auto &cell : tria.active_cell_iterators())
+      {
+        switch (cell->active_cell_index())
+          {
+            case 0:
+              cell->set_subdomain_id(0);
+              break;
+            case 1:
+            case 2:
+              cell->set_subdomain_id(2);
+              break;
+            case 3:
+              cell->set_subdomain_id(3);
+              break;
+          }
+      }
+  });
+
+  GridGenerator::subdivided_hyper_cube(tria, 2);
+
+  const auto construction_data =
+    TriangulationDescription::Utilities::create_description_from_triangulation(
+      tria,
+      RepartitioningPolicyTools::DefaultPolicy<dim>(true).partition(tria));
+
+  parallel::fullydistributed::Triangulation<dim> tria_pft(comm);
+  tria_pft.create_triangulation(construction_data);
+
+  deallog << tria.n_locally_owned_active_cells() << " "
+          << tria_pft.n_locally_owned_active_cells() << std::endl;
+}
+
+
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  MPI_Comm comm = MPI_COMM_WORLD;
+
+  deallog.push("all");
+  test<2>(comm);
+  deallog.pop();
+}
diff --git a/tests/fullydistributed_grids/repartitioning_07.with_p4est=true.mpirun=4.output b/tests/fullydistributed_grids/repartitioning_07.with_p4est=true.mpirun=4.output
new file mode 100644 (file)
index 0000000..17d8de7
--- /dev/null
@@ -0,0 +1,11 @@
+
+DEAL:0:all::1 1
+
+DEAL:1:all::0 2
+
+
+DEAL:2:all::2 1
+
+
+DEAL:3:all::1 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.