]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add function 8444/head
authortcclevenger <tcleven@clemson.edu>
Sat, 3 Aug 2019 01:06:36 +0000 (19:06 -0600)
committertcclevenger <tcleven@clemson.edu>
Tue, 6 Aug 2019 15:12:04 +0000 (09:12 -0600)
include/deal.II/multigrid/mg_tools.h
source/multigrid/mg_tools.cc
source/multigrid/mg_tools.inst.in
tests/multigrid/workload_imbalance.cc [new file with mode: 0644]
tests/multigrid/workload_imbalance.with_mpi=true.with_p4est=true.mpirun=3.output [new file with mode: 0644]

index 99048d1ecf04685e5ce4d51c402ec0af0d9956c4..c38bfaf58b79e9da7ac58e22ac9ee7b0901e6fa0 100644 (file)
@@ -279,6 +279,19 @@ namespace MGTools
   template <int dim, int spacedim>
   unsigned int
   max_level_for_coarse_mesh(const Triangulation<dim, spacedim> &tria);
+
+  /**
+   * Return the imbalance of the parallel distribution of the multigrid
+   * mesh hierarchy. Ideally this value is equal to 1 (every processor owns
+   * the same number of cells on each level, approximately true for most
+   * globally refined meshes). Values greater than 1 estimate the slowdown
+   * one should see in a geometric multigrid v-cycle as compared with the same
+   * computation on a perfectly distributed mesh hierarchy.
+   */
+  template <int dim, int spacedim>
+  double
+  workload_imbalance(const Triangulation<dim, spacedim> &tria);
+
 } // namespace MGTools
 
 /* @} */
index d107f40bda0a59ccff1d00fc10dbb4169baf7332..4d3c8c04f2571142a8bdb9f97147c3a99fe534e9 100644 (file)
@@ -1661,6 +1661,58 @@ namespace MGTools
 
     return global_min;
   }
+
+
+
+  template <int dim, int spacedim>
+  double
+  workload_imbalance(const Triangulation<dim, spacedim> &tria)
+  {
+    double workload_imbalance = 1.0;
+
+    // It is only necessary to calculate the imbalance
+    // on a distributed mesh. The imbalance is always
+    // 1.0 for the serial case.
+    if (const parallel::Triangulation<dim, spacedim> *tr =
+          dynamic_cast<const parallel::Triangulation<dim, spacedim> *>(&tria))
+      {
+        const unsigned int n_proc =
+          Utilities::MPI::n_mpi_processes(tr->get_communicator());
+        const unsigned int n_global_levels = tr->n_global_levels();
+
+        // This value will represent the sum over the multigrid
+        // levels of the maximum number of cells owned by any
+        // one processesor on that level.
+        types::global_dof_index work_estimate = 0;
+
+        // Sum of all cells in the multigrid hierarchy
+        types::global_dof_index total_cells_in_hierarchy = 0;
+
+        for (int lvl = n_global_levels - 1; lvl >= 0; --lvl)
+          {
+            // Number of cells this processor owns on this level
+            types::global_dof_index n_owned_cells_on_lvl = 0;
+
+            for (const auto &cell : tr->cell_iterators_on_level(lvl))
+              if (cell->is_locally_owned_on_level())
+                ++n_owned_cells_on_lvl;
+
+            work_estimate +=
+              dealii::Utilities::MPI::max(n_owned_cells_on_lvl,
+                                          tr->get_communicator());
+
+            total_cells_in_hierarchy +=
+              dealii::Utilities::MPI::sum(n_owned_cells_on_lvl,
+                                          tr->get_communicator());
+          }
+
+        const double ideal_work =
+          total_cells_in_hierarchy / static_cast<double>(n_proc);
+        workload_imbalance = work_estimate / ideal_work;
+      }
+
+    return workload_imbalance;
+  }
 } // namespace MGTools
 
 
index 1a78b8e2ef2a2af514d30a995af697d5d4053d6b..e1252193fb3228a23cfb7c0aaed1829569a3f1ee 100644 (file)
@@ -162,6 +162,10 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
       template unsigned int
       max_level_for_coarse_mesh(
         const Triangulation<deal_II_dimension, deal_II_space_dimension> &tria);
+
+      template double
+      workload_imbalance(
+        const Triangulation<deal_II_dimension, deal_II_space_dimension> &tria);
 #endif
     \}
   }
diff --git a/tests/multigrid/workload_imbalance.cc b/tests/multigrid/workload_imbalance.cc
new file mode 100644 (file)
index 0000000..20492c9
--- /dev/null
@@ -0,0 +1,81 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2006 - 2018 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.
+//
+// ---------------------------------------------------------------------
+
+
+// check function MGTools::workload_imbalance()
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/distributed/grid_refinement.h>
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+
+#include <deal.II/multigrid/mg_constrained_dofs.h>
+#include <deal.II/multigrid/multigrid.h>
+
+#include <deal.II/numerics/vector_tools.h>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+test()
+{
+  parallel::distributed::Triangulation<dim> tria(
+    MPI_COMM_WORLD,
+    Triangulation<dim>::limit_level_difference_at_vertices,
+    parallel::distributed::Triangulation<dim>::construct_multigrid_hierarchy);
+  GridGenerator::hyper_cube(tria, 0, 1);
+  tria.refine_global(1);
+
+  unsigned int counter = 0;
+  for (auto cell : tria.active_cell_iterators())
+    {
+      cell->set_refine_flag();
+      ++counter;
+
+      if (dim == 2 && counter == 1)
+        break;
+      if (dim == 3 && counter == 4)
+        break;
+    }
+  tria.execute_coarsening_and_refinement();
+
+  const double workload_imbalance = MGTools::workload_imbalance(tria);
+  deallog << "Workload imbalance: " << workload_imbalance << std::endl;
+}
+
+
+
+int
+main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  deallog.push("2d");
+  test<2>();
+  deallog.pop();
+  deallog.push("3d");
+  test<3>();
+  deallog.pop();
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/multigrid/workload_imbalance.with_mpi=true.with_p4est=true.mpirun=3.output b/tests/multigrid/workload_imbalance.with_mpi=true.with_p4est=true.mpirun=3.output
new file mode 100644 (file)
index 0000000..9d5b5fc
--- /dev/null
@@ -0,0 +1,14 @@
+
+DEAL:0:2d::Workload imbalance: 2.66667
+DEAL:0:3d::Workload imbalance: 1.60976
+DEAL:0::OK
+
+DEAL:1:2d::Workload imbalance: 2.66667
+DEAL:1:3d::Workload imbalance: 1.60976
+DEAL:1::OK
+
+
+DEAL:2:2d::Workload imbalance: 2.66667
+DEAL:2:3d::Workload imbalance: 1.60976
+DEAL:2::OK
+

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.