]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use Utilties::MPI::compute_index_owner() instead of lower-level facilities. 17704/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 13 Sep 2024 13:27:53 +0000 (07:27 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Fri, 13 Sep 2024 13:35:42 +0000 (07:35 -0600)
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
source/distributed/repartitioning_policy_tools.cc
source/multigrid/mg_tools.cc

index 7db8a1b7f1981ed62f6d788f06b7e88f5728753a..eca23769822af62a663f13de4acfa8b258f0eff2 100644 (file)
@@ -837,25 +837,11 @@ namespace internal
 
           is_dst_remote_potentially_relevant.subtract_set(is_dst_locally_owned);
 
-          std::vector<unsigned int> owning_ranks_of_ghosts(
-            is_dst_remote_potentially_relevant.n_elements());
-
-          {
-            Utilities::MPI::internal::ComputeIndexOwner::
-              ConsensusAlgorithmsPayload process(
-                is_dst_locally_owned,
-                is_dst_remote_potentially_relevant,
-                communicator,
-                owning_ranks_of_ghosts,
-                false);
-
-            Utilities::MPI::ConsensusAlgorithms::Selector<
-              std::vector<
-                std::pair<types::global_cell_index, types::global_cell_index>>,
-              std::vector<unsigned int>>
-              consensus_algorithm;
-            consensus_algorithm.run(process, communicator);
-          }
+          const std::vector<unsigned int> owning_ranks_of_ghosts =
+            Utilities::MPI::compute_index_owner(
+              is_dst_locally_owned,
+              is_dst_remote_potentially_relevant,
+              communicator);
 
           for (unsigned i = 0;
                i < is_dst_remote_potentially_relevant.n_elements();
@@ -3961,22 +3947,10 @@ MGTwoLevelTransfer<dim, VectorType>::reinit(
 
       const MPI_Comm communicator = dof_handler_fine.get_mpi_communicator();
 
-      std::vector<unsigned int> owning_ranks(
-        is_locally_owned_coarse.n_elements());
-
-      Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload
-        process(is_locally_owned_fine,
-                is_locally_owned_coarse,
-                communicator,
-                owning_ranks,
-                false);
-
-      Utilities::MPI::ConsensusAlgorithms::Selector<
-        std::vector<
-          std::pair<types::global_cell_index, types::global_cell_index>>,
-        std::vector<unsigned int>>
-        consensus_algorithm;
-      consensus_algorithm.run(process, communicator);
+      const std::vector<unsigned int> owning_ranks =
+        Utilities::MPI::compute_index_owner(is_locally_owned_fine,
+                                            is_locally_owned_coarse,
+                                            communicator);
 
       bool all_cells_found = true;
 
index 1377c7d1a45f82703c069862b1ec0feb4eb197d0..7430a8fc36f0a9ce02d6032ab8b02845e2d0e156 100644 (file)
@@ -13,8 +13,7 @@
 // ------------------------------------------------------------------------
 
 
-#include <deal.II/base/mpi_compute_index_owner_internal.h>
-#include <deal.II/base/mpi_consensus_algorithms.h>
+#include <deal.II/base/mpi.h>
 
 #include <deal.II/distributed/repartitioning_policy_tools.h>
 #include <deal.II/distributed/tria_base.h>
@@ -145,28 +144,14 @@ namespace RepartitioningPolicyTools
       if (cell->is_locally_owned())
         is_coarse.add_index(cell_id_translator.translate(cell));
 
-    std::vector<unsigned int> owning_ranks_of_coarse_cells(
-      is_coarse.n_elements());
-    {
-      Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload
-        process(is_level_partitions,
-                is_coarse,
-                communicator,
-                owning_ranks_of_coarse_cells,
-                false);
-
-      Utilities::MPI::ConsensusAlgorithms::Selector<
-        std::vector<
-          std::pair<types::global_cell_index, types::global_cell_index>>,
-        std::vector<unsigned int>>
-        consensus_algorithm;
-      consensus_algorithm.run(process, communicator);
-    }
+    const std::vector<unsigned int> owning_ranks_of_coarse_cells =
+      Utilities::MPI::compute_index_owner(is_level_partitions,
+                                          is_coarse,
+                                          communicator);
 
     const auto tria =
       dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
         &tria_coarse_in);
-
     Assert(tria, ExcNotImplemented());
 
     LinearAlgebra::distributed::Vector<double> partition(
index 62ecfa943e699a486190d15802a88c13c746c3e3..db27943f7450b67723c8cedfe26484a6602714d8 100644 (file)
@@ -14,7 +14,7 @@
 
 #include <deal.II/base/logstream.h>
 #include <deal.II/base/mg_level_object.h>
-#include <deal.II/base/mpi_compute_index_owner_internal.h>
+#include <deal.II/base/mpi.h>
 #include <deal.II/base/thread_management.h>
 
 #include <deal.II/distributed/tria_base.h>
@@ -1771,22 +1771,10 @@ namespace MGTools
                   cell_id_translator.translate(cell, i));
             }
 
-        std::vector<unsigned int> is_fine_required_ranks(
-          is_fine_required.n_elements());
-
-        Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload
-          process(is_fine_owned,
-                  is_fine_required,
-                  communicator,
-                  is_fine_required_ranks,
-                  false);
-
-        Utilities::MPI::ConsensusAlgorithms::Selector<
-          std::vector<
-            std::pair<types::global_cell_index, types::global_cell_index>>,
-          std::vector<unsigned int>>
-          consensus_algorithm;
-        consensus_algorithm.run(process, communicator);
+        const std::vector<unsigned int> is_fine_required_ranks =
+          Utilities::MPI::compute_index_owner(is_fine_owned,
+                                              is_fine_required,
+                                              communicator);
 
         for (unsigned i = 0; i < is_fine_required.n_elements(); ++i)
           if (is_fine_required_ranks[i] == my_rank)

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.