]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add default constructor to NoncontiguousPartitioner 9633/head
authorPeter Munch <peterrmuench@gmail.com>
Sat, 7 Mar 2020 10:01:05 +0000 (11:01 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Sun, 8 Mar 2020 06:21:40 +0000 (07:21 +0100)
include/deal.II/base/mpi_noncontiguous_partitioner.h
include/deal.II/base/mpi_noncontiguous_partitioner.templates.h

index e04ac51190de8f8d9ae4c9d2b6c3a629482f2baa..1942348c49b263fcfe1af6193657acf5a63dc34f 100644 (file)
@@ -43,6 +43,12 @@ namespace Utilities
       : public LinearAlgebra::CommunicationPatternBase
     {
     public:
+      /**
+       * Default constructor. Requires calling one of the reinit() functions
+       * to create a valid object.
+       */
+      NoncontiguousPartitioner() = default;
+
       /**
        * Constructor. Set up point-to-point communication pattern based on the
        * IndexSets arguments @p indexset_has and @p indexset_want for the MPI
@@ -124,6 +130,14 @@ namespace Utilities
              const IndexSet &indexset_want,
              const MPI_Comm &communicator) override;
 
+      /**
+       * Initialize the inner data structures.
+       */
+      void
+      reinit(const std::vector<types::global_dof_index> &indices_has,
+             const std::vector<types::global_dof_index> &indices_want,
+             const MPI_Comm &                            communicator);
+
     private:
       /// MPI communicator
       MPI_Comm communicator;
index e8b6ca134ff40bd92167bb7095c5fb121218acb7..9fe98cc0dbb00d41a67ad0eb26830a583bac9d8a 100644 (file)
@@ -42,60 +42,18 @@ namespace Utilities
       this->reinit(indexset_has, indexset_want, communicator);
     }
 
+
+
     template <typename Number>
     NoncontiguousPartitioner<Number>::NoncontiguousPartitioner(
       const std::vector<types::global_dof_index> &indices_has,
       const std::vector<types::global_dof_index> &indices_want,
       const MPI_Comm &                            communicator)
     {
-      // step 0) determine "number of degrees of freedom" needed for IndexSet
-      const types::global_dof_index n_dofs = Utilities::MPI::max(
-        std::max(
-          [&indices_has]() {
-            const auto it = max_element(indices_has.begin(), indices_has.end());
-            return it == indices_has.end() ? 0 : (*it + 1);
-          }(),
-          [&indices_want]() {
-            const auto it =
-              max_element(indices_want.begin(), indices_want.end());
-            return it == indices_want.end() ? 0 : (*it + 1);
-          }()),
-        communicator);
-
-      // step 1) convert vectors to indexsets (sorted!)
-      IndexSet index_set_has(n_dofs);
-      index_set_has.add_indices(indices_has.begin(), indices_has.end());
-
-      IndexSet index_set_want(n_dofs);
-      index_set_want.add_indices(indices_want.begin(), indices_want.end());
-
-      // step 2) setup internal data structures with indexset
-      this->reinit(index_set_has, index_set_want, communicator);
-
-      // step 3) fix inner data structures so that it is sorted as
-      //         in the original vector
-      {
-        std::vector<types::global_dof_index> temp_map_send(
-          index_set_has.n_elements());
-
-        for (types::global_dof_index i = 0; i < indices_has.size(); i++)
-          temp_map_send[index_set_has.index_within_set(indices_has[i])] = i;
-
-        for (auto &i : send_indices)
-          i = temp_map_send[i];
-      }
-
-      {
-        std::vector<types::global_dof_index> temp_map_recv(
-          index_set_want.n_elements());
+      this->reinit(indices_has, indices_want, communicator);
+    }
 
-        for (types::global_dof_index i = 0; i < indices_want.size(); i++)
-          temp_map_recv[index_set_want.index_within_set(indices_want[i])] = i;
 
-        for (auto &i : recv_indices)
-          i = temp_map_recv[i];
-      }
-    }
 
     template <typename Number>
     std::pair<unsigned int, unsigned int>
@@ -104,6 +62,8 @@ namespace Utilities
       return {send_ranks.size(), recv_ranks.size()};
     }
 
+
+
     template <typename Number>
     types::global_dof_index
     NoncontiguousPartitioner<Number>::memory_consumption()
@@ -120,6 +80,8 @@ namespace Utilities
              MemoryConsumption::memory_consumption(recv_requests);
     }
 
+
+
     template <typename Number>
     const MPI_Comm &
     NoncontiguousPartitioner<Number>::get_mpi_communicator() const
@@ -127,6 +89,8 @@ namespace Utilities
       return communicator;
     }
 
+
+
     template <typename Number>
     void
     NoncontiguousPartitioner<Number>::reinit(const IndexSet &indexset_has,
@@ -212,6 +176,63 @@ namespace Utilities
 
 
 
+    template <typename Number>
+    void
+    NoncontiguousPartitioner<Number>::reinit(
+      const std::vector<types::global_dof_index> &indices_has,
+      const std::vector<types::global_dof_index> &indices_want,
+      const MPI_Comm &                            communicator)
+    {
+      // step 0) determine "number of degrees of freedom" needed for IndexSet
+      const types::global_dof_index n_dofs = Utilities::MPI::max(
+        std::max(
+          [&indices_has]() {
+            const auto it = max_element(indices_has.begin(), indices_has.end());
+            return it == indices_has.end() ? 0 : (*it + 1);
+          }(),
+          [&indices_want]() {
+            const auto it =
+              max_element(indices_want.begin(), indices_want.end());
+            return it == indices_want.end() ? 0 : (*it + 1);
+          }()),
+        communicator);
+
+      // step 1) convert vectors to indexsets (sorted!)
+      IndexSet index_set_has(n_dofs);
+      index_set_has.add_indices(indices_has.begin(), indices_has.end());
+
+      IndexSet index_set_want(n_dofs);
+      index_set_want.add_indices(indices_want.begin(), indices_want.end());
+
+      // step 2) setup internal data structures with indexset
+      this->reinit(index_set_has, index_set_want, communicator);
+
+      // step 3) fix inner data structures so that it is sorted as
+      //         in the original vector
+      {
+        std::vector<types::global_dof_index> temp_map_send(
+          index_set_has.n_elements());
+
+        for (types::global_dof_index i = 0; i < indices_has.size(); i++)
+          temp_map_send[index_set_has.index_within_set(indices_has[i])] = i;
+
+        for (auto &i : send_indices)
+          i = temp_map_send[i];
+      }
+
+      {
+        std::vector<types::global_dof_index> temp_map_recv(
+          index_set_want.n_elements());
+
+        for (types::global_dof_index i = 0; i < indices_want.size(); i++)
+          temp_map_recv[index_set_want.index_within_set(indices_want[i])] = i;
+
+        for (auto &i : recv_indices)
+          i = temp_map_recv[i];
+      }
+    }
+
+
     template <typename Number>
     template <typename VectorType>
     void

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.