]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Avoid double packing/unpacking for CA algorithms. 13733/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 13 May 2022 22:31:04 +0000 (16:31 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Fri, 13 May 2022 22:31:04 +0000 (16:31 -0600)
include/deal.II/fe/fe_tools_extrapolate.templates.h
source/grid/grid_tools.cc
source/grid/tria_description.cc
source/multigrid/mg_tools.cc

index c883d07ad2143d4f4aeb07ccb81563137e530ffc..fc92f2a2a6eca38dff00b9646bbf5d658e767b8b 100644 (file)
@@ -20,6 +20,7 @@
 #include <deal.II/base/config.h>
 
 #include <deal.II/base/mpi_consensus_algorithms.h>
+#include <deal.II/base/mpi_consensus_algorithms.templates.h>
 
 #include <deal.II/distributed/p4est_wrappers.h>
 #include <deal.II/distributed/tria.h>
@@ -1140,45 +1141,45 @@ namespace FETools
       // front which processes (and from how many processes) we have to
       // expect information from.
       const auto create_request =
-        [&cells_to_send](const types::subdomain_id other_rank) {
-          std::vector<CellData> cells_for_this_destination;
-          for (const auto &cell : cells_to_send)
-            if (cell.receiver == other_rank)
-              cells_for_this_destination.emplace_back(cell);
-
-          return Utilities::pack(cells_for_this_destination, false);
-        };
+        [&cells_to_send](
+          const types::subdomain_id other_rank) -> std::vector<CellData> {
+        std::vector<CellData> cells_for_this_destination;
+        for (const auto &cell : cells_to_send)
+          if (cell.receiver == other_rank)
+            cells_for_this_destination.emplace_back(cell);
+
+        return cells_for_this_destination;
+      };
 
       const auto answer_request =
-        [&received_cells](const unsigned int       other_rank,
-                          const std::vector<char> &request) {
-          // We got a message from 'other_rank', so let us decode the
-          // message in the same way as we have assembled it above.
-          // Note that the cells just received do not contain
-          // information where they came from, and we have to add that
-          // ourselves for later use.
-          for (CellData &cell_data :
-               Utilities::unpack<std::vector<CellData>>(request, false))
-            {
-              cell_data.receiver = other_rank;
-              received_cells.emplace_back(std::move(cell_data));
-            }
+        [&received_cells](const unsigned int           other_rank,
+                          const std::vector<CellData> &request) -> int {
+        // We got a message from 'other_rank', so let us decode the
+        // message in the same way as we have assembled it above.
+        // Note that the cells just received do not contain
+        // information where they came from, and we have to add that
+        // ourselves for later use.
+        for (CellData cell_data : request)
+          {
+            cell_data.receiver = other_rank;
+            received_cells.emplace_back(std::move(cell_data));
+          }
 
-          // Nothing left to do here, we don't actually need to provide an
-          // answer:
-          return std::vector<char>();
-        };
+        // Nothing left to do here, we don't actually need to provide an
+        // answer:
+        return 0;
+      };
 
       const auto read_answer = [](const unsigned int /*other_rank*/,
-                                  const std::vector<char> &answer) {
+                                  const int &answer) {
         // We don't put anything into the answers, so nothing should
-        // have been coming out at this end either:
+        // have been coming out at this end either that differs from
+        // the default-sent zero integer:
         (void)answer;
-        Assert(answer.size() == 0, ExcInternalError());
+        Assert(answer == 0, ExcInternalError());
       };
 
-      Utilities::MPI::ConsensusAlgorithms::selector<std::vector<char>,
-                                                    std::vector<char>>(
+      Utilities::MPI::ConsensusAlgorithms::selector<std::vector<CellData>, int>(
         destinations,
         create_request,
         answer_request,
index 1ddaf4030a69327f0e523840cd8f8054442f74d1..85ba1aa885af184efa846bf510b8c6a3b757d799 100644 (file)
@@ -16,6 +16,7 @@
 #include <deal.II/base/mpi.h>
 #include <deal.II/base/mpi.templates.h>
 #include <deal.II/base/mpi_consensus_algorithms.h>
+#include <deal.II/base/mpi_consensus_algorithms.templates.h>
 #include <deal.II/base/quadrature_lib.h>
 #include <deal.II/base/thread_management.h>
 
@@ -58,6 +59,7 @@
 
 #include <deal.II/physics/transformations.h>
 
+
 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
 #include <boost/random/mersenne_twister.hpp>
 #include <boost/random/uniform_real_distribution.hpp>
@@ -5985,6 +5987,9 @@ namespace GridTools
           "The marked_vertices vector has to be either empty or its size has "
           "to equal the number of vertices of the triangulation."));
 
+      using RequestType = std::vector<std::pair<unsigned int, Point<spacedim>>>;
+      using AnswerType  = std::vector<unsigned int>;
+
       // In the case that a marked_vertices vector has been given and none
       // of its entries is true, we know that this process does not own
       // any of the incoming points (and it will not send any data) so
@@ -5997,36 +6002,31 @@ namespace GridTools
       const auto create_request = [&](const unsigned int other_rank) {
         const auto other_rank_index = translate(other_rank);
 
-        std::vector<std::pair<unsigned int, Point<spacedim>>> temp;
-        temp.reserve(potential_owners_ptrs[other_rank_index + 1] -
-                     potential_owners_ptrs[other_rank_index]);
+        RequestType request;
+        request.reserve(potential_owners_ptrs[other_rank_index + 1] -
+                        potential_owners_ptrs[other_rank_index]);
 
         for (unsigned int i = potential_owners_ptrs[other_rank_index];
              i < potential_owners_ptrs[other_rank_index + 1];
              ++i)
-          temp.emplace_back(potential_owners_indices[i],
-                            points[potential_owners_indices[i]]);
+          request.emplace_back(potential_owners_indices[i],
+                               points[potential_owners_indices[i]]);
 
-        return Utilities::pack(temp, false);
+        return request;
       };
 
       const auto answer_request =
-        [&](const unsigned int &     other_rank,
-            const std::vector<char> &request) -> std::vector<char> {
-        const auto recv_buffer_unpacked = Utilities::unpack<
-          std::vector<std::pair<unsigned int, Point<spacedim>>>>(request,
-                                                                 false);
-
-        std::vector<unsigned int> request_buffer_temp(
-          recv_buffer_unpacked.size(), 0);
+        [&](const unsigned int &other_rank,
+            const RequestType & request) -> AnswerType {
+        AnswerType answer(request.size(), 0);
 
         if (has_relevant_vertices)
           {
             cell_hint = cache.get_triangulation().begin_active();
 
-            for (unsigned int i = 0; i < recv_buffer_unpacked.size(); ++i)
+            for (unsigned int i = 0; i < request.size(); ++i)
               {
-                const auto &index_and_point = recv_buffer_unpacked[i];
+                const auto &index_and_point = request[i];
 
                 const auto cells_and_reference_positions =
                   find_all_locally_owned_active_cells_around_point(
@@ -6051,27 +6051,24 @@ namespace GridTools
                       numbers::invalid_unsigned_int);
                   }
 
-                request_buffer_temp[i] = cells_and_reference_positions.size();
+                answer[i] = cells_and_reference_positions.size();
               }
           }
 
         if (perform_handshake)
-          return Utilities::pack(request_buffer_temp, false);
+          return answer;
         else
           return {};
       };
 
-      const auto process_answer = [&](const unsigned int       other_rank,
-                                      const std::vector<char> &answer) {
+      const auto process_answer = [&](const unsigned int other_rank,
+                                      const AnswerType & answer) {
         if (perform_handshake)
           {
-            const auto recv_buffer_unpacked =
-              Utilities::unpack<std::vector<unsigned int>>(answer, false);
-
             const auto other_rank_index = translate(other_rank);
 
-            for (unsigned int i = 0; i < recv_buffer_unpacked.size(); ++i)
-              for (unsigned int j = 0; j < recv_buffer_unpacked[i]; ++j)
+            for (unsigned int i = 0; i < answer.size(); ++i)
+              for (unsigned int j = 0; j < answer[i]; ++j)
                 recv_components.emplace_back(
                   other_rank,
                   potential_owners_indices
@@ -6080,8 +6077,7 @@ namespace GridTools
           }
       };
 
-      Utilities::MPI::ConsensusAlgorithms::selector<std::vector<char>,
-                                                    std::vector<char>>(
+      Utilities::MPI::ConsensusAlgorithms::selector<RequestType, AnswerType>(
         potential_owners_ranks,
         create_request,
         answer_request,
index aea46c970a2f055675053a7e99e4f7c658b2797a..c58ef335ebec394ef50209fe82a1804b802f0f64 100644 (file)
@@ -16,6 +16,7 @@
 #include <deal.II/base/geometry_info.h>
 #include <deal.II/base/mpi.h>
 #include <deal.II/base/mpi_consensus_algorithms.h>
+#include <deal.II/base/mpi_consensus_algorithms.templates.h>
 
 #include <deal.II/distributed/fully_distributed_tria.h>
 #include <deal.II/distributed/tria.h>
@@ -123,32 +124,19 @@ namespace TriangulationDescription
             const auto other_rank_index =
               std::distance(relevant_processes.begin(), ptr);
 
-            return dealii::Utilities::pack(description_temp[other_rank_index],
-                                           false);
-          };
-
-          const auto answer_request = [&](const unsigned int,
-                                          const std::vector<char> &request) {
-            this->merge(
-              dealii::Utilities::unpack<DescriptionTemp<dim, spacedim>>(request,
-                                                                        false),
-              vertices_have_unique_ids);
-            return std::vector<char>();
-          };
-
-          const auto process_answer = [](const unsigned int,
-                                         const std::vector<char> &answer) {
-            (void)answer;
-            Assert(answer.size() == 0, ExcInternalError());
+            return description_temp[other_rank_index];
           };
 
+          const auto answer_request =
+            [&](const unsigned int,
+                const DescriptionTemp<dim, spacedim> &request) {
+              this->merge(request, vertices_have_unique_ids);
+              return char{};
+            };
 
-          dealii::Utilities::MPI::ConsensusAlgorithms::
-            selector<std::vector<char>, std::vector<char>>(relevant_processes,
-                                                           create_request,
-                                                           answer_request,
-                                                           process_answer,
-                                                           comm);
+          dealii::Utilities::MPI::ConsensusAlgorithms::selector<
+            DescriptionTemp<dim, spacedim>,
+            char>(relevant_processes, create_request, answer_request, {}, comm);
         }
 
         /**
index e1a9c19405cf1be3ae4fed41a91d6cbd8c870d3a..b139e56a1d129666e2efd04374447a0155f368b4 100644 (file)
@@ -16,6 +16,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_consensus_algorithms.templates.h>
 #include <deal.II/base/thread_management.h>
 
 #include <deal.II/dofs/dof_accessor.h>

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.