]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Compute number of locally owned active cells per process on-demand. 8808/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 20 Sep 2019 06:43:52 +0000 (08:43 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Fri, 20 Sep 2019 20:56:12 +0000 (22:56 +0200)
17 files changed:
doc/news/changes/incompatibilities/20190920PeterMunch [new file with mode: 0644]
include/deal.II/distributed/tria_base.h
source/distributed/tria_base.cc
tests/mpi/cell_weights_01.cc
tests/mpi/cell_weights_01_back_and_forth_01.cc
tests/mpi/cell_weights_01_back_and_forth_02.cc
tests/mpi/cell_weights_02.cc
tests/mpi/cell_weights_03.cc
tests/mpi/cell_weights_04.cc
tests/mpi/cell_weights_05.cc
tests/mpi/cell_weights_06.cc
tests/mpi/hp_step-40.cc
tests/mpi/hp_step-40_variable_01.cc
tests/mpi/step-40.cc
tests/mpi/step-40_cuthill_mckee.cc
tests/mpi/step-40_cuthill_mckee_MPI-subset.cc
tests/mpi/step-40_direct_solver.cc

diff --git a/doc/news/changes/incompatibilities/20190920PeterMunch b/doc/news/changes/incompatibilities/20190920PeterMunch
new file mode 100644 (file)
index 0000000..d830fdf
--- /dev/null
@@ -0,0 +1,5 @@
+Changed: The parallel::TriangulationBase class does not store the number of active cells
+of all MPI processes any more. Instead, the information is computed on-demand when calling
+the function parallel::TriangulationBase::compute_n_locally_owned_active_cells_per_processor.
+<br>
+(Peter Munch, 2019/09/20)
index 53749beb9938448a4a0e7c4d48eb8cfe8b35efc9..fafcd5d1c2429cce26e6c3ba13716b4c96d92d29 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2008 - 2018 by the deal.II authors
+// Copyright (C) 2008 - 2019 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -79,10 +79,11 @@ namespace parallel
      * that contribute to this triangulation. The element of this vector
      * indexed by locally_owned_subdomain() equals the result of
      * n_locally_owned_active_cells().
+     *
+     * @note This function involves global communication!
      */
-    const std::vector<unsigned int> &
-    n_locally_owned_active_cells_per_processor() const;
-
+    std::vector<unsigned int>
+    compute_n_locally_owned_active_cells_per_processor() const;
 
     /**
      * Return the number of active cells in the triangulation that are locally
@@ -195,10 +196,9 @@ namespace parallel
     struct NumberCache
     {
       /**
-       * This vector stores the number of locally owned active cells per MPI
-       * rank.
+       * Number of locally owned active cells of this MPI rank.
        */
-      std::vector<unsigned int> n_locally_owned_active_cells;
+      unsigned int n_locally_owned_active_cells;
       /**
        * The total number of active cells (sum of @p
        * n_locally_owned_active_cells).
index 44252fa12f287f8aefa03f70cf322b3525c84e37..b3c6d4b5c7ceb3accadbc48961898d7bd28ddf1d 100644 (file)
@@ -56,7 +56,6 @@ namespace parallel
            ExcMessage("You compiled deal.II without MPI support, for "
                       "which parallel::TriangulationBase is not available."));
 #endif
-    number_cache.n_locally_owned_active_cells.resize(n_subdomains);
   }
 
 
@@ -98,8 +97,6 @@ namespace parallel
       this->dealii::Triangulation<dim, spacedim>::memory_consumption() +
       MemoryConsumption::memory_consumption(mpi_communicator) +
       MemoryConsumption::memory_consumption(my_subdomain) +
-      MemoryConsumption::memory_consumption(
-        number_cache.n_locally_owned_active_cells) +
       MemoryConsumption::memory_consumption(
         number_cache.n_global_active_cells) +
       MemoryConsumption::memory_consumption(number_cache.n_global_levels);
@@ -117,7 +114,7 @@ namespace parallel
 
   template <int dim, int spacedim>
   TriangulationBase<dim, spacedim>::NumberCache::NumberCache()
-    : n_global_active_cells(0)
+    : n_locally_owned_active_cells(0)
     , n_global_levels(0)
   {}
 
@@ -125,7 +122,7 @@ namespace parallel
   unsigned int
   TriangulationBase<dim, spacedim>::n_locally_owned_active_cells() const
   {
-    return number_cache.n_locally_owned_active_cells[my_subdomain];
+    return number_cache.n_locally_owned_active_cells;
   }
 
   template <int dim, int spacedim>
@@ -143,11 +140,32 @@ namespace parallel
   }
 
   template <int dim, int spacedim>
-  const std::vector<unsigned int> &
-  TriangulationBase<dim, spacedim>::n_locally_owned_active_cells_per_processor()
-    const
+  std::vector<unsigned int>
+  TriangulationBase<dim, spacedim>::
+    compute_n_locally_owned_active_cells_per_processor() const
   {
-    return number_cache.n_locally_owned_active_cells;
+    ;
+#ifdef DEAL_II_WITH_MPI
+    std::vector<unsigned int> n_locally_owned_active_cells_per_processor(
+      Utilities::MPI::n_mpi_processes(this->mpi_communicator), 0);
+
+    if (this->n_levels() > 0)
+      {
+        const int ierr =
+          MPI_Allgather(&number_cache.n_locally_owned_active_cells,
+                        1,
+                        MPI_UNSIGNED,
+                        n_locally_owned_active_cells_per_processor.data(),
+                        1,
+                        MPI_UNSIGNED,
+                        this->mpi_communicator);
+        AssertThrowMPI(ierr);
+      }
+
+    return n_locally_owned_active_cells_per_processor;
+#else
+    return {number_cache.n_locally_owned_active_cells};
+#endif
   }
 
   template <int dim, int spacedim>
@@ -162,16 +180,9 @@ namespace parallel
   void
   TriangulationBase<dim, spacedim>::update_number_cache()
   {
-    Assert(number_cache.n_locally_owned_active_cells.size() ==
-             Utilities::MPI::n_mpi_processes(this->mpi_communicator),
-           ExcInternalError());
-
-    std::fill(number_cache.n_locally_owned_active_cells.begin(),
-              number_cache.n_locally_owned_active_cells.end(),
-              0);
-
     number_cache.ghost_owners.clear();
     number_cache.level_ghost_owners.clear();
+    number_cache.n_locally_owned_active_cells = 0;
 
     if (this->n_levels() == 0)
       {
@@ -206,25 +217,11 @@ namespace parallel
            cell != this->end();
            ++cell)
         if (cell->subdomain_id() == my_subdomain)
-          ++number_cache.n_locally_owned_active_cells[my_subdomain];
-
-    unsigned int send_value =
-      number_cache.n_locally_owned_active_cells[my_subdomain];
-    const int ierr =
-      MPI_Allgather(&send_value,
-                    1,
-                    MPI_UNSIGNED,
-                    number_cache.n_locally_owned_active_cells.data(),
-                    1,
-                    MPI_UNSIGNED,
-                    this->mpi_communicator);
-    AssertThrowMPI(ierr);
+          ++number_cache.n_locally_owned_active_cells;
 
     number_cache.n_global_active_cells =
-      std::accumulate(number_cache.n_locally_owned_active_cells.begin(),
-                      number_cache.n_locally_owned_active_cells.end(),
-                      /* ensure sum is computed with correct data type:*/
-                      static_cast<types::global_dof_index>(0));
+      Utilities::MPI::sum(number_cache.n_locally_owned_active_cells,
+                          this->mpi_communicator);
     number_cache.n_global_levels =
       Utilities::MPI::max(this->n_levels(), this->mpi_communicator);
   }
index 32fa0bac5e80853688fb94f22c966ee0ab8f19e7..932e9a28588af6cd49e91c3bc45bf9c22f148843 100644 (file)
@@ -53,10 +53,12 @@ test()
   // (roughly) equal between all processors
   tr.repartition();
 
+  const auto n_locally_owned_active_cells_per_processor =
+    tr.compute_n_locally_owned_active_cells_per_processor();
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     for (unsigned int p = 0; p < numproc; ++p)
       deallog << "processor " << p << ": "
-              << tr.n_locally_owned_active_cells_per_processor()[p]
+              << n_locally_owned_active_cells_per_processor[p]
               << " locally owned active cells" << std::endl;
 }
 
index 9ebbc79689304b10186ca546ae85b45d91300029..2686080140c18754b3387a3742a1c4b7244aaa9e 100644 (file)
@@ -87,11 +87,12 @@ test()
                                            std::placeholders::_2));
   tr.repartition();
 
-
+  const auto n_locally_owned_active_cells_per_processor =
+    tr.compute_n_locally_owned_active_cells_per_processor();
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     for (unsigned int p = 0; p < numproc; ++p)
       deallog << "processor " << p << ": "
-              << tr.n_locally_owned_active_cells_per_processor()[p]
+              << n_locally_owned_active_cells_per_processor[p]
               << " locally owned active cells" << std::endl;
 }
 
index 4d28e2a0b2aada603c37cfd9e31686bde264933d..edbffd902f23f4ea38f1c275298356c5dd0a7c39 100644 (file)
@@ -72,11 +72,12 @@ test()
   tr.signals.cell_weight.disconnect_all_slots();
   tr.repartition();
 
-
+  const auto n_locally_owned_active_cells_per_processor =
+    tr.compute_n_locally_owned_active_cells_per_processor();
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     for (unsigned int p = 0; p < numproc; ++p)
       deallog << "processor " << p << ": "
-              << tr.n_locally_owned_active_cells_per_processor()[p]
+              << n_locally_owned_active_cells_per_processor[p]
               << " locally owned active cells" << std::endl;
 }
 
index 20fd4cd4a357f49937484757a386f7e651f694f6..862b0644fc8ed8353985b035a7996c43d816da26 100644 (file)
@@ -59,10 +59,12 @@ test()
     std::bind(&cell_weight<dim>, std::placeholders::_1, std::placeholders::_2));
   tr.refine_global(1);
 
+  const auto n_locally_owned_active_cells_per_processor =
+    tr.compute_n_locally_owned_active_cells_per_processor();
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     for (unsigned int p = 0; p < numproc; ++p)
       deallog << "processor " << p << ": "
-              << tr.n_locally_owned_active_cells_per_processor()[p]
+              << n_locally_owned_active_cells_per_processor[p]
               << " locally owned active cells" << std::endl;
 }
 
index de63dea8e36a9b4d30292ce50133239734c7c9d0..5a44ca85eba0ae0c3dae8c4673947b61d455f145 100644 (file)
@@ -67,10 +67,12 @@ test()
 
   tr.repartition();
 
+  const auto n_locally_owned_active_cells_per_processor =
+    tr.compute_n_locally_owned_active_cells_per_processor();
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     for (unsigned int p = 0; p < numproc; ++p)
       deallog << "processor " << p << ": "
-              << tr.n_locally_owned_active_cells_per_processor()[p]
+              << n_locally_owned_active_cells_per_processor[p]
               << " locally owned active cells" << std::endl;
 
   // let each processor sum up its weights
index 367eb75caff988a3cb39c33c7ef7115231f5a3c9..fab4b09040b4e5482b53140bce7715066ffa295c 100644 (file)
@@ -61,10 +61,12 @@ test()
 
   tr.refine_global(1);
 
+  const auto n_locally_owned_active_cells_per_processor =
+    tr.compute_n_locally_owned_active_cells_per_processor();
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     for (unsigned int p = 0; p < numproc; ++p)
       deallog << "processor " << p << ": "
-              << tr.n_locally_owned_active_cells_per_processor()[p]
+              << n_locally_owned_active_cells_per_processor[p]
               << " locally owned active cells" << std::endl;
 
   // let each processor sum up its weights
index df14695434b44eac7c97be7287d4be3b9366c137..b0c2bee4d9c5d05bcac462c7396beaa20d8d9331 100644 (file)
@@ -80,10 +80,12 @@ test()
 
   tr.repartition();
 
+  const auto n_locally_owned_active_cells_per_processor =
+    tr.compute_n_locally_owned_active_cells_per_processor();
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     for (unsigned int p = 0; p < numproc; ++p)
       deallog << "processor " << p << ": "
-              << tr.n_locally_owned_active_cells_per_processor()[p]
+              << n_locally_owned_active_cells_per_processor[p]
               << " locally owned active cells" << std::endl;
 
   // let each processor sum up its weights
index 3d1e6b7f4a2dead1d4f2d8bbdc2d04e98757fab1..b57d4f3f102fd74a0d3e33643dd77d1824267f6d 100644 (file)
@@ -79,10 +79,12 @@ test()
     std::bind(&cell_weight<dim>, std::placeholders::_1, std::placeholders::_2));
   tr.repartition();
 
+  const auto n_locally_owned_active_cells_per_processor =
+    tr.compute_n_locally_owned_active_cells_per_processor();
   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
     for (unsigned int p = 0; p < numproc; ++p)
       deallog << "processor " << p << ": "
-              << tr.n_locally_owned_active_cells_per_processor()[p]
+              << n_locally_owned_active_cells_per_processor[p]
               << " locally owned active cells" << std::endl;
 
   // let each processor sum up its weights
index ba44608f2cc2d5c4f8bc67b58990b92d79274a49..44918b2fb0ae60d79cebbe4c872f6e0bb0953a1b 100644 (file)
@@ -329,11 +329,12 @@ namespace Step40
         pcout << "   Number of active cells:       "
               << triangulation.n_global_active_cells() << std::endl
               << "      ";
+        const auto n_locally_owned_active_cells_per_processor =
+          triangulation.compute_n_locally_owned_active_cells_per_processor();
         for (unsigned int i = 0;
              i < Utilities::MPI::n_mpi_processes(mpi_communicator);
              ++i)
-          pcout << triangulation.n_locally_owned_active_cells_per_processor()[i]
-                << '+';
+          pcout << n_locally_owned_active_cells_per_processor[i] << '+';
         pcout << std::endl;
 
         pcout << "   Number of degrees of freedom: " << dof_handler.n_dofs()
index 73e117f835b5cd20a19164922266b966b009a593..90427a37a73f1f344b54ab0edfc0d6a8a3405059 100644 (file)
@@ -332,11 +332,12 @@ namespace Step40
         pcout << "   Number of active cells:       "
               << triangulation.n_global_active_cells() << std::endl
               << "      ";
+        const auto n_locally_owned_active_cells_per_processor =
+          triangulation.compute_n_locally_owned_active_cells_per_processor();
         for (unsigned int i = 0;
              i < Utilities::MPI::n_mpi_processes(mpi_communicator);
              ++i)
-          pcout << triangulation.n_locally_owned_active_cells_per_processor()[i]
-                << '+';
+          pcout << n_locally_owned_active_cells_per_processor[i] << '+';
         pcout << std::endl;
 
         pcout << "   Number of degrees of freedom: " << dof_handler.n_dofs()
index 6253d4b7a54d0306e6deee1b28f610726ef4091c..092bb521b830b5b77a3b7a623385864ca5f29551 100644 (file)
@@ -314,11 +314,12 @@ namespace Step40
         pcout << "   Number of active cells:       "
               << triangulation.n_global_active_cells() << std::endl
               << "      ";
+        const auto n_locally_owned_active_cells_per_processor =
+          triangulation.compute_n_locally_owned_active_cells_per_processor();
         for (unsigned int i = 0;
              i < Utilities::MPI::n_mpi_processes(mpi_communicator);
              ++i)
-          pcout << triangulation.n_locally_owned_active_cells_per_processor()[i]
-                << '+';
+          pcout << n_locally_owned_active_cells_per_processor[i] << '+';
         pcout << std::endl;
 
         pcout << "   Number of degrees of freedom: " << dof_handler.n_dofs()
index 5225ccaaa3acc9698d9e9e15a930d45c7c5539c4..65a0289b108a2e3f615a870ce67ad79f8ea2a08b 100644 (file)
@@ -378,11 +378,12 @@ namespace Step40
         pcout << "   Number of active cells:       "
               << triangulation.n_global_active_cells() << std::endl
               << "      ";
+        const auto n_locally_owned_active_cells_per_processor =
+          triangulation.compute_n_locally_owned_active_cells_per_processor();
         for (unsigned int i = 0;
              i < Utilities::MPI::n_mpi_processes(mpi_communicator);
              ++i)
-          pcout << triangulation.n_locally_owned_active_cells_per_processor()[i]
-                << '+';
+          pcout << n_locally_owned_active_cells_per_processor[i] << '+';
         pcout << std::endl;
 
         pcout << "   Number of degrees of freedom: " << dof_handler.n_dofs()
index b6582551261dc613a38f232e681919b87605ff21..7516e42d8066417e6b9f5436a3ab7f08ec1d4fcb 100644 (file)
@@ -379,11 +379,12 @@ namespace Step40
         pcout << "   Number of active cells:       "
               << triangulation.n_global_active_cells() << std::endl
               << "      ";
+        const auto n_locally_owned_active_cells_per_processor =
+          triangulation.compute_n_locally_owned_active_cells_per_processor();
         for (unsigned int i = 0;
              i < Utilities::MPI::n_mpi_processes(mpi_communicator);
              ++i)
-          pcout << triangulation.n_locally_owned_active_cells_per_processor()[i]
-                << '+';
+          pcout << n_locally_owned_active_cells_per_processor[i] << '+';
         pcout << std::endl;
 
         pcout << "   Number of degrees of freedom: " << dof_handler.n_dofs()
index 6b80fca5c7802c34c4766726646c65a34829f464..1a43bb5cee4361755c462d03a3b6cc03c506e499 100644 (file)
@@ -292,11 +292,12 @@ namespace Step40
         pcout << "   Number of active cells:       "
               << triangulation.n_global_active_cells() << std::endl
               << "      ";
+        const auto n_locally_owned_active_cells_per_processor =
+          triangulation.compute_n_locally_owned_active_cells_per_processor();
         for (unsigned int i = 0;
              i < Utilities::MPI::n_mpi_processes(mpi_communicator);
              ++i)
-          pcout << triangulation.n_locally_owned_active_cells_per_processor()[i]
-                << '+';
+          pcout << n_locally_owned_active_cells_per_processor[i] << '+';
         pcout << std::endl;
 
         pcout << "   Number of degrees of freedom: " << dof_handler.n_dofs()

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.