]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove template argument from LinearAlgebra::SharedMPI::Partitioner
authorPeter Munch <peterrmuench@gmail.com>
Mon, 23 Mar 2020 09:47:56 +0000 (10:47 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 23 Mar 2020 09:47:56 +0000 (10:47 +0100)
include/deal.II/lac/la_sm_partitioner.h
include/deal.II/lac/la_sm_vector.h
include/deal.II/lac/la_sm_vector.templates.h
include/deal.II/matrix_free/matrix_free.h

index 2b60c5e4fd541196ddc649f36e02267c74b30dfb..7973662c901fdb6be0e20aa9b3ed372cd81a9151 100644 (file)
@@ -30,7 +30,6 @@ namespace LinearAlgebra
 {
   namespace SharedMPI
   {
-    template <typename Number = double>
     class Partitioner : public LinearAlgebra::CommunicationPatternBase
     {
     public:
@@ -228,18 +227,21 @@ namespace LinearAlgebra
                       send_sm_req.data(),
                       MPI_STATUSES_IGNORE);
         }
-
-        buffer.resize(send_remote_ptr.back());
       }
 
+      template <typename Number>
       void
       update_ghost_values_start(
-        Number *               data_this,
-        std::vector<Number *> &data_others,
-        const unsigned int     communication_channel = 0) const
+        Number *                       data_this,
+        std::vector<Number *> &        data_others,
+        dealii::AlignedVector<Number> &buffer,
+        const unsigned int             communication_channel = 0) const
       {
         (void)data_others;
 
+        if (send_remote_ptr.back() != buffer.size())
+          buffer.resize(send_remote_ptr.back());
+
         int dummy;
         for (unsigned int i = 0; i < recv_sm_ranks.size(); i++)
           MPI_Irecv(&dummy,
@@ -285,10 +287,14 @@ namespace LinearAlgebra
           }
       }
 
+      template <typename Number>
       void
-      update_ghost_values_finish(Number *               data_this,
-                                 std::vector<Number *> &data_others) const
+      update_ghost_values_finish(Number *                       data_this,
+                                 std::vector<Number *> &        data_others,
+                                 dealii::AlignedVector<Number> &buffer) const
       {
+        (void)buffer;
+
         for (unsigned int c = 0; c < recv_sm_ranks.size(); c++)
           {
             int i;
@@ -318,21 +324,28 @@ namespace LinearAlgebra
                     MPI_STATUSES_IGNORE);
       }
 
+      template <typename Number>
       void
-      update_ghost_values(Number *               data_this,
-                          std::vector<Number *> &data_others) const
+      update_ghost_values(Number *                       data_this,
+                          std::vector<Number *> &        data_others,
+                          dealii::AlignedVector<Number> &buffer) const
       {
-        update_ghost_values_start(data_this, data_others);
-        update_ghost_values_finish(data_this, data_others);
+        update_ghost_values_start(data_this, data_others, buffer);
+        update_ghost_values_finish(data_this, data_others, buffer);
       }
 
+      template <typename Number>
       void
-      compress_start(Number *               data_this,
-                     std::vector<Number *> &data_others,
-                     const unsigned int     communication_channel = 0) const
+      compress_start(Number *                       data_this,
+                     std::vector<Number *> &        data_others,
+                     dealii::AlignedVector<Number> &buffer,
+                     const unsigned int communication_channel = 0) const
       {
         (void)data_others;
 
+        if (send_remote_ptr.back() != buffer.size())
+          buffer.resize(send_remote_ptr.back());
+
         int dummy;
         for (unsigned int i = 0; i < recv_sm_ranks.size(); i++)
           MPI_Isend(&dummy,
@@ -371,9 +384,11 @@ namespace LinearAlgebra
                     send_remote_req.data() + i);
       }
 
+      template <typename Number>
       void
-      compress_finish(Number *               data_this,
-                      std::vector<Number *> &data_others) const
+      compress_finish(Number *                       data_this,
+                      std::vector<Number *> &        data_others,
+                      dealii::AlignedVector<Number> &buffer) const
       {
         for (unsigned int c = 0; c < send_sm_ranks.size(); c++)
           {
@@ -419,11 +434,14 @@ namespace LinearAlgebra
                     MPI_STATUSES_IGNORE);
       }
 
+      template <typename Number>
       void
-      compress(Number *data_this, std::vector<Number *> &data_others) const
+      compress(Number *                       data_this,
+               std::vector<Number *> &        data_others,
+               dealii::AlignedVector<Number> &buffer) const
       {
-        compress_start(data_this, data_others);
-        compress_finish(data_this, data_others);
+        compress_start(data_this, data_others, buffer);
+        compress_finish(data_this, data_others, buffer);
       }
 
     private:
@@ -432,8 +450,6 @@ namespace LinearAlgebra
 
       unsigned int n_local_elements;
 
-      mutable AlignedVector<Number> buffer;
-
       std::vector<unsigned int>            recv_remote_ranks;
       std::vector<types::global_dof_index> recv_remote_ptr = {0};
       mutable std::vector<MPI_Request>     recv_remote_req;
index 4e89827b0162b9ed98edda6293d771a46f1a68ac..423e3403e181591c92764867fb149cf0a615a4a3 100644 (file)
@@ -159,7 +159,7 @@ namespace LinearAlgebra
       void
       reinit(
         const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
-        const std::shared_ptr<const Partitioner<Number>> &partitioner_sm);
+        const std::shared_ptr<const Partitioner> &partitioner_sm);
 
       void
       swap(Vector<Number, MemorySpace> &v);
@@ -477,7 +477,7 @@ namespace LinearAlgebra
 
       std::shared_ptr<const Utilities::MPI::Partitioner> partitioner;
 
-      std::shared_ptr<const Partitioner<Number>> partitioner_sm;
+      std::shared_ptr<const Partitioner> partitioner_sm;
 
       size_type allocated_size;
 
@@ -491,7 +491,7 @@ namespace LinearAlgebra
         thread_loop_partitioner;
 
       // needed?
-      mutable MemorySpaceData<Number> import_data;
+      mutable dealii::AlignedVector<Number> import_data;
 
       mutable bool vector_is_ghosted;
 
index 41841b283a93ecdae1e14a63e07abb6b841a1aea..796cfff466c73c73f84d2e61700b4167cc468d90 100644 (file)
@@ -206,8 +206,7 @@ namespace LinearAlgebra
       // is only used as temporary storage for compress() and
       // update_ghost_values, and we might have vectors where we never
       // call these methods and hence do not need to have the storage.
-      import_data.values.reset();
-      import_data.values_dev.reset();
+      import_data.clear(); // TODO
 
       thread_loop_partitioner = v.thread_loop_partitioner;
     }
@@ -246,7 +245,7 @@ namespace LinearAlgebra
     void
     Vector<Number, MemorySpaceType>::reinit(
       const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
-      const std::shared_ptr<const Partitioner<Number>> &        partitioner_sm)
+      const std::shared_ptr<const Partitioner> &                partitioner_sm)
     {
       clear_mpi_requests();
       this->partitioner    = partitioner;
@@ -265,8 +264,7 @@ namespace LinearAlgebra
       // is only used as temporary storage for compress() and
       // update_ghost_values, and we might have vectors where we never
       // call these methods and hence do not need to have the storage.
-      import_data.values.reset();
-      import_data.values_dev.reset();
+      import_data.clear();
 
       vector_is_ghosted = false;
     }
@@ -523,6 +521,7 @@ namespace LinearAlgebra
              ExcMessage("Cannot call compress() on a ghosted vector"));
       partitioner_sm->compress_start(data.values.get(),
                                      data.others,
+                                     import_data,
                                      communication_channel);
     }
 
@@ -536,7 +535,9 @@ namespace LinearAlgebra
       Assert(::dealii::VectorOperation::values::add == operation,
              ExcNotImplemented());
       vector_is_ghosted = false;
-      partitioner_sm->compress_finish(data.values.get(), data.others);
+      partitioner_sm->compress_finish(data.values.get(),
+                                      data.others,
+                                      import_data);
     }
 
 
@@ -548,6 +549,7 @@ namespace LinearAlgebra
     {
       partitioner_sm->update_ghost_values_start(data.values.get(),
                                                 data.others,
+                                                import_data,
                                                 communication_channel);
     }
 
@@ -558,7 +560,8 @@ namespace LinearAlgebra
     Vector<Number, MemorySpaceType>::update_ghost_values_finish() const
     {
       partitioner_sm->update_ghost_values_finish(data.values.get(),
-                                                 data.others);
+                                                 data.others,
+                                                 import_data);
       vector_is_ghosted = true;
     }
 
index b4509d9c6772dd4eee8cfb06c0237755909efcd9..922bb61b1dfd2ab6fc573f42ac39ed54441f4836 100644 (file)
@@ -2095,9 +2095,8 @@ private:
   // TODO: move it to DoFInfo
   mutable std::map<
     unsigned int,
-    std::map<
-      const MPI_Comm,
-      std::shared_ptr<const LinearAlgebra::SharedMPI::Partitioner<Number>>>>
+    std::map<const MPI_Comm,
+             std::shared_ptr<const LinearAlgebra::SharedMPI::Partitioner>>>
     partitioner_sm;
 
   /**
@@ -2240,7 +2239,7 @@ MatrixFree<dim, Number, VectorizedArrayType>::initialize_dof_vector(
 
   if (partitioner_sm[comp][comm_sm] == nullptr)
     partitioner_sm[comp][comm_sm] =
-      std::make_shared<LinearAlgebra::SharedMPI::Partitioner<Number>>(
+      std::make_shared<LinearAlgebra::SharedMPI::Partitioner>(
         part->get_communicator(),
         comm_sm,
         part->locally_owned_range(),

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.