]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add empty partitioner
authorPeter Munch <peterrmuench@gmail.com>
Sat, 21 Mar 2020 13:32:41 +0000 (14:32 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Sat, 21 Mar 2020 13:32:41 +0000 (14:32 +0100)
include/deal.II/lac/la_sm_partitioner.h [new file with mode: 0644]
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
tests/sm/vector_sm_01.cc

diff --git a/include/deal.II/lac/la_sm_partitioner.h b/include/deal.II/lac/la_sm_partitioner.h
new file mode 100644 (file)
index 0000000..5245b31
--- /dev/null
@@ -0,0 +1,67 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_la_sm_partitioner_h
+#define dealii_la_sm_partitioner_h
+
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/lac/communication_pattern_base.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+
+namespace LinearAlgebra
+{
+  namespace SharedMPI
+  {
+    template <typename Number = double>
+    class Partitioner : public LinearAlgebra::CommunicationPatternBase
+    {
+    public:
+      Partitioner(const MPI_Comm &comm, const MPI_Comm &comm_sm)
+        : comm(comm)
+        , comm_sm(comm_sm)
+      {}
+
+      const MPI_Comm &
+      get_mpi_communicator() const override
+      {
+        return comm;
+      }
+
+      void
+      reinit(const IndexSet &indexset_has,
+             const IndexSet &indexset_want,
+             const MPI_Comm &communicator) override
+      {
+        (void)indexset_has;
+        (void)indexset_want;
+        (void)communicator;
+      }
+
+    private:
+      const MPI_Comm &comm;
+      const MPI_Comm &comm_sm;
+    };
+
+  } // end of namespace SharedMPI
+} // end of namespace LinearAlgebra
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
\ No newline at end of file
index 748b73db318cee8a60a5158b3658bdcdb96548a1..5a60b1dfab060950cc2ef64c98ca258fbeeb2dab 100644 (file)
@@ -24,6 +24,7 @@
 #include <deal.II/base/partitioner.h>
 #include <deal.II/base/thread_management.h>
 
+#include <deal.II/lac/la_sm_partitioner.h>
 #include <deal.II/lac/vector_operation.h>
 #include <deal.II/lac/vector_space_vector.h>
 #include <deal.II/lac/vector_type_traits.h>
@@ -157,7 +158,8 @@ namespace LinearAlgebra
 
       void
       reinit(
-        const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner);
+        const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
+        const std::shared_ptr<const Partitioner<Number>> &partitioner_sm);
 
       void
       swap(Vector<Number, MemorySpace> &v);
@@ -475,6 +477,8 @@ namespace LinearAlgebra
 
       std::shared_ptr<const Utilities::MPI::Partitioner> partitioner;
 
+      std::shared_ptr<const Partitioner<Number>> partitioner_sm;
+
       size_type allocated_size;
 
       mutable MemorySpaceData<Number> data;
index fd4ca5de3b43844a77158395603eca5640ca5324..2e46004f13f75b71b38fd51596c1001ce0ffab52 100644 (file)
@@ -13,8 +13,8 @@
 //
 // ---------------------------------------------------------------------
 
-#ifndef dealii_la_parallel_vector_templates_h
-#define dealii_la_parallel_vector_templates_h
+#ifndef dealii_la_sm_vector_templates_h
+#define dealii_la_sm_vector_templates_h
 
 
 #include <deal.II/base/config.h>
@@ -170,14 +170,16 @@ namespace LinearAlgebra
     template <typename Number, typename MemorySpaceType>
     void
     Vector<Number, MemorySpaceType>::reinit(
-      const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_in)
+      const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
+      const std::shared_ptr<const Partitioner<Number>> &        partitioner_sm)
     {
       clear_mpi_requests();
-      partitioner = partitioner_in;
+      this->partitioner    = partitioner;
+      this->partitioner_sm = partitioner_sm;
 
       // set vector size and allocate memory
       const size_type new_allocated_size =
-        partitioner->local_size() + partitioner->n_ghost_indices();
+        this->partitioner->local_size() + this->partitioner->n_ghost_indices();
       resize_val(new_allocated_size);
 
       // initialize to zero
index 082d5492137469de4c85fdbdbd5abc522a92f5e2..a3d9078ea1d6c29f1ea73eb8cefeafeaaae6ba42 100644 (file)
@@ -40,6 +40,7 @@
 #include <deal.II/lac/affine_constraints.h>
 #include <deal.II/lac/block_vector_base.h>
 #include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/la_sm_partitioner.h>
 #include <deal.II/lac/la_sm_vector.h>
 #include <deal.II/lac/vector_operation.h>
 
@@ -1481,6 +1482,7 @@ public:
   template <typename Number2>
   void
   initialize_dof_vector(LinearAlgebra::SharedMPI::Vector<Number2> &vec,
+                        const MPI_Comm                             comm_sm,
                         const unsigned int dof_handler_index = 0) const;
 
   /**
@@ -2090,6 +2092,14 @@ private:
    */
   std::vector<internal::MatrixFreeFunctions::DoFInfo> dof_info;
 
+  // TODO: move it to DoFInfo
+  mutable std::map<
+    unsigned int,
+    std::map<
+      const MPI_Comm,
+      std::shared_ptr<const LinearAlgebra::SharedMPI::Partitioner<Number>>>>
+    partitioner_sm;
+
   /**
    * Contains the weights for constraints stored in DoFInfo. Filled into a
    * separate field since several vector components might share similar
@@ -2221,10 +2231,17 @@ template <typename Number2>
 inline void
 MatrixFree<dim, Number, VectorizedArrayType>::initialize_dof_vector(
   LinearAlgebra::SharedMPI::Vector<Number2> &vec,
+  const MPI_Comm                             comm_sm,
   const unsigned int                         comp) const
 {
   AssertIndexRange(comp, n_components());
-  vec.reinit(dof_info[comp].vector_partitioner);
+
+  if (partitioner_sm[comp][comm_sm] == nullptr)
+    partitioner_sm[comp][comm_sm] =
+      std::make_shared<LinearAlgebra::SharedMPI::Partitioner<Number>>(
+        dof_info[comp].vector_partitioner->get_communicator(), comm_sm);
+
+  vec.reinit(dof_info[comp].vector_partitioner, partitioner_sm[comp][comm_sm]);
 }
 
 
index 749c8d033d387de3f5d1d30b39cb8dfc8080a998..17d609739178c9320469042a233de4afc51018a2 100644 (file)
@@ -60,8 +60,10 @@ test(const int n_refinements, const int degree)
   dealii::MatrixFree<dim, Number> matrix_free;
   matrix_free.reinit(mapping, dof_handler, constraint, quad, additional_data);
 
+  MPI_Comm comm_sm;
+
   LinearAlgebra::SharedMPI::Vector<Number> vec;
-  matrix_free.initialize_dof_vector(vec);
+  matrix_free.initialize_dof_vector(vec, comm_sm);
 }
 
 int

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.