]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add Triangulation/DoFHandler::get_communicator() 11822/head
authorPeter Munch <peterrmuench@gmail.com>
Sun, 28 Feb 2021 09:39:45 +0000 (10:39 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 1 Mar 2021 16:55:41 +0000 (17:55 +0100)
doc/news/changes/major/20210228Munch [new file with mode: 0644]
include/deal.II/distributed/tria_base.h
include/deal.II/dofs/dof_handler.h
include/deal.II/grid/tria.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
source/distributed/tria_base.cc
source/grid/tria.cc

diff --git a/doc/news/changes/major/20210228Munch b/doc/news/changes/major/20210228Munch
new file mode 100644 (file)
index 0000000..dc73da8
--- /dev/null
@@ -0,0 +1,6 @@
+New: The communicator of an arbitrary (not just parallel) Triangulation class can now be 
+queried via Triangulation::get_communicator() or DoFHandler::get_communicator(). In
+the case of serial Triangulations and DoFHandler set up with serial Triangulations,
+MPI_COMM_SELF is returned.
+<br>
+(Peter Munch, 2021/02/28)
index 19bcf9f8dc596083a7b007c16595d93aa8685183..e11f58764fa19421514bb93fab26f2eb6710b1f7 100644 (file)
@@ -96,8 +96,8 @@ namespace parallel
     /**
      * Return MPI communicator used by this triangulation.
      */
-    virtual const MPI_Comm &
-    get_communicator() const;
+    virtual MPI_Comm
+    get_communicator() const override;
 
     /**
      * Return if multilevel hierarchy is supported and has been constructed.
@@ -308,7 +308,7 @@ namespace parallel
      * communicator for this class, which is a duplicate of the one passed to
      * the constructor.
      */
-    MPI_Comm mpi_communicator;
+    const MPI_Comm mpi_communicator;
 
     /**
      * The subdomain id to be used for the current processor. This is the MPI
index c74286fd7f0d37e1acab01f722671ff32ca80e4a..99601c45c7e136182d4552b151d7e58f66af22d1 100644 (file)
@@ -1267,6 +1267,12 @@ public:
   const Triangulation<dim, spacedim> &
   get_triangulation() const;
 
+  /**
+   * Return MPI communicator used by the underlying triangulation.
+   */
+  MPI_Comm
+  get_communicator() const;
+
   /**
    * Whenever serialization with a parallel::distributed::Triangulation as the
    * underlying triangulation is considered, we also need to consider storing
@@ -2019,6 +2025,18 @@ DoFHandler<dim, spacedim>::get_triangulation() const
 
 
 
+template <int dim, int spacedim>
+inline MPI_Comm
+DoFHandler<dim, spacedim>::get_communicator() const
+{
+  Assert(tria != nullptr,
+         ExcMessage("This DoFHandler object has not been associated "
+                    "with a triangulation."));
+  return tria->get_communicator();
+}
+
+
+
 template <int dim, int spacedim>
 inline const BlockInfo &
 DoFHandler<dim, spacedim>::block_info() const
index c8e4eb1166c89c1ae0222b7c12593718535d0e37..a840993ab3bf0941b67bdce238c0be1efb97d9c5 100644 (file)
@@ -1614,6 +1614,13 @@ public:
   virtual void
   clear();
 
+  /**
+   * Return MPI communicator used by this triangulation. In the case of
+   * a serial Triangulation object, MPI_COMM_SELF is returned.
+   */
+  virtual MPI_Comm
+  get_communicator() const;
+
   /**
    * Set the mesh smoothing to @p mesh_smoothing. This overrides the
    * MeshSmoothing given to the constructor. It is allowed to call this
index 933f550cbc508f10cfc9ca380a9ced8cee010dd7..43e9856705fc103d98594c82eaa15f44939ed8f1 100644 (file)
@@ -238,19 +238,6 @@ namespace internal
 {
   namespace
   {
-    template <typename MeshType>
-    MPI_Comm
-    get_mpi_comm(const MeshType &mesh)
-    {
-      const auto *tria_parallel = dynamic_cast<
-        const dealii::parallel::TriangulationBase<MeshType::dimension,
-                                                  MeshType::space_dimension> *>(
-        &(mesh.get_triangulation()));
-
-      return tria_parallel != nullptr ? tria_parallel->get_communicator() :
-                                        MPI_COMM_SELF;
-    }
-
     template <int dim>
     unsigned int
     compute_shift_within_children(const unsigned int child,
@@ -516,7 +503,8 @@ namespace internal
     FineDoFHandlerView(const MeshType &mesh_fine, const MeshType &mesh_coarse)
       : mesh_fine(mesh_fine)
       , mesh_coarse(mesh_coarse)
-      , communicator(get_mpi_comm(mesh_fine) /*TODO: fix for different comms*/)
+      , communicator(
+          mesh_fine.get_communicator() /*TODO: fix for different comms*/)
       , cell_id_translator(n_coarse_cells(mesh_fine),
                            n_global_levels(mesh_fine))
     {
@@ -873,7 +861,7 @@ namespace internal
           n_coarse_cells =
             std::max(n_coarse_cells, cell->id().get_coarse_cell_id());
 
-      return Utilities::MPI::max(n_coarse_cells, get_mpi_comm(mesh)) + 1;
+      return Utilities::MPI::max(n_coarse_cells, mesh.get_communicator()) + 1;
     }
 
     static unsigned int
@@ -973,9 +961,10 @@ namespace internal
               std::max(max_active_fe_indices[1], cell->active_fe_index());
           }
 
-      const auto comm = get_mpi_comm(dof_handler_fine);
+      const auto comm = dof_handler_fine.get_communicator();
 
-      Assert(comm == get_mpi_comm(dof_handler_coarse), ExcNotImplemented());
+      Assert(comm == dof_handler_coarse.get_communicator(),
+             ExcNotImplemented());
 
       ArrayView<unsigned int> temp_min(min_active_fe_indices);
       ArrayView<unsigned int> temp_max(max_active_fe_indices);
@@ -1002,10 +991,10 @@ namespace internal
       {
         // ... for fine mesh
         {
-          transfer.partitioner_fine.reset(
-            new Utilities::MPI::Partitioner(view.locally_owned_dofs(),
-                                            view.locally_relevant_dofs(),
-                                            get_mpi_comm(dof_handler_fine)));
+          transfer.partitioner_fine.reset(new Utilities::MPI::Partitioner(
+            view.locally_owned_dofs(),
+            view.locally_relevant_dofs(),
+            dof_handler_fine.get_communicator()));
           transfer.vec_fine.reinit(transfer.partitioner_fine);
         }
 
@@ -1018,7 +1007,7 @@ namespace internal
           transfer.partitioner_coarse.reset(new Utilities::MPI::Partitioner(
             dof_handler_coarse.locally_owned_dofs(),
             locally_relevant_dofs,
-            get_mpi_comm(dof_handler_coarse)));
+            dof_handler_coarse.get_communicator()));
           transfer.vec_coarse.reinit(transfer.partitioner_coarse);
         }
       }
@@ -1324,7 +1313,7 @@ namespace internal
             std::make_shared<Utilities::MPI::Partitioner>(
               dof_handler_fine.locally_owned_dofs(),
               locally_relevant_dofs,
-              get_mpi_comm(dof_handler_fine));
+              dof_handler_fine.get_communicator());
           transfer.vec_fine.reinit(transfer.partitioner_fine);
 
           touch_count_.reinit(partitioner_fine_);
@@ -1487,7 +1476,7 @@ namespace internal
 
       transfer.constraint_coarse.copy_from(constraint_coarse);
 
-      const auto comm = get_mpi_comm(dof_handler_coarse);
+      const auto comm = dof_handler_coarse.get_communicator();
       {
         IndexSet locally_relevant_dofs;
         DoFTools::extract_locally_relevant_dofs(dof_handler_fine,
index a27a626cee34db7b96a6e530a882b9d3c9e5bde8..b379bb8ffa27e5a6c68fcf4b685682c2a7793d3c 100644 (file)
@@ -90,7 +90,7 @@ namespace parallel
   {
     std::size_t mem =
       this->dealii::Triangulation<dim, spacedim>::memory_consumption() +
-      MemoryConsumption::memory_consumption(mpi_communicator) +
+      MemoryConsumption::memory_consumption(this->mpi_communicator) +
       MemoryConsumption::memory_consumption(my_subdomain) +
       MemoryConsumption::memory_consumption(
         number_cache.n_global_active_cells) +
@@ -135,7 +135,7 @@ namespace parallel
   }
 
   template <int dim, int spacedim>
-  const MPI_Comm &
+  MPI_Comm
   TriangulationBase<dim, spacedim>::get_communicator() const
   {
     return mpi_communicator;
@@ -170,7 +170,7 @@ namespace parallel
           number_cache.ghost_owners.insert(cell->subdomain_id());
 
       Assert(number_cache.ghost_owners.size() <
-               Utilities::MPI::n_mpi_processes(mpi_communicator),
+               Utilities::MPI::n_mpi_processes(this->mpi_communicator),
              ExcInternalError());
     }
 
index 1fd8e850c127d533313ad28f6764e50cb6bda31f..7884f80b19cb959f7d805aa7b2110da753b105ae 100644 (file)
@@ -10190,6 +10190,14 @@ Triangulation<dim, spacedim>::clear()
 }
 
 
+template <int dim, int spacedim>
+MPI_Comm
+Triangulation<dim, spacedim>::get_communicator() const
+{
+  return MPI_COMM_SELF;
+}
+
+
 
 template <int dim, int spacedim>
 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.