]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce DistributedTrinagulationBase 8592/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 16 Aug 2019 05:27:23 +0000 (07:27 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Thu, 22 Aug 2019 03:46:10 +0000 (05:46 +0200)
Conflicts:
include/deal.II/distributed/tria.h
include/deal.II/distributed/tria_base.h
source/distributed/tria.cc

doc/news/changes/minor/20190817PeterMunch [new file with mode: 0644]
include/deal.II/distributed/tria.h
include/deal.II/distributed/tria_base.h
include/deal.II/dofs/dof_handler_policy.h
source/distributed/tria.cc
source/distributed/tria_base.cc
source/distributed/tria_base.inst.in
source/dofs/dof_handler_policy.cc

diff --git a/doc/news/changes/minor/20190817PeterMunch b/doc/news/changes/minor/20190817PeterMunch
new file mode 100644 (file)
index 0000000..80e03e2
--- /dev/null
@@ -0,0 +1,4 @@
+New: Introduce new class parallel::DistributedTriangulationBase 
+between parallel::TriangulationBase and parallel::distributed::Triangulation.
+<br>
+(Peter Munch, 2019/08/17)
index cb70df8a589d08c57889c8f39fd6d0190509a0df..3737ca094aa1c0cf1a19aff3cd569d9959e9bb4f 100644 (file)
@@ -249,7 +249,7 @@ namespace parallel
      */
     template <int dim, int spacedim = dim>
     class Triangulation
-      : public dealii::parallel::TriangulationBase<dim, spacedim>
+      : public dealii::parallel::DistributedTriangulationBase<dim, spacedim>
     {
     public:
       /**
@@ -383,6 +383,12 @@ namespace parallel
       virtual void
       clear() override;
 
+      /**
+       * Return if multilevel hierarchy is supported and has been constructed.
+       */
+      bool
+      is_multilevel_hierarchy_constructed() const override;
+
       /**
        * Implementation of the same function as in the base class.
        *
@@ -1247,7 +1253,7 @@ namespace parallel
      */
     template <int spacedim>
     class Triangulation<1, spacedim>
-      : public dealii::parallel::TriangulationBase<1, spacedim>
+      : public dealii::parallel::DistributedTriangulationBase<1, spacedim>
     {
     public:
       /**
@@ -1335,6 +1341,9 @@ namespace parallel
       void
       save(const std::string &filename) const;
 
+      bool
+      is_multilevel_hierarchy_constructed() const override;
+
       /**
        * This function is not implemented, but needs to be present for the
        * compiler.
index 09f5fa3e3ac0f8e00104679dd2b1ed26c9c7866d..53749beb9938448a4a0e7c4d48eb8cfe8b35efc9 100644 (file)
@@ -246,6 +246,32 @@ namespace parallel
   template <int dim, int spacedim = dim>
   using Triangulation DEAL_II_DEPRECATED = TriangulationBase<dim, spacedim>;
 
+  /**
+   * A base class for distributed triangulations. It can be used to test
+   * (via <tt>dynamic_cast</tt>) whether a particular triangulation object
+   * is able to distribute meshes between processes.
+   */
+  template <int dim, int spacedim = dim>
+  class DistributedTriangulationBase
+    : public dealii::parallel::TriangulationBase<dim, spacedim>
+  {
+  public:
+    /**
+     * Constructor.
+     */
+    DistributedTriangulationBase(
+      MPI_Comm mpi_communicator,
+      const typename dealii::Triangulation<dim, spacedim>::MeshSmoothing
+                 smooth_grid = (dealii::Triangulation<dim, spacedim>::none),
+      const bool check_for_distorted_cells = false);
+
+    /**
+     * Return if multilevel hierarchy is supported and has been constructed.
+     */
+    virtual bool
+    is_multilevel_hierarchy_constructed() const = 0;
+  };
+
 } // namespace parallel
 
 DEAL_II_NAMESPACE_CLOSE
index 9f08546f6128e13f133d86c4935c7f3aa4d80e08..4c6616705014865d5af78c8735523b3b6c2c86ef 100644 (file)
@@ -215,7 +215,7 @@ namespace internal
 
       /**
        * This class implements the policy for operations when we use a
-       * parallel::distributed::Triangulation object.
+       * parallel::DistributedTriangulationBase object.
        */
       template <class DoFHandlerType>
       class ParallelDistributed
index a186866e3df7da7ccdc15c79bce319967b1e5888..58d00ebd913277f58944141b9b61117bf1c5988a 100644 (file)
@@ -2127,7 +2127,7 @@ namespace parallel
       : // Do not check for distorted cells.
         // For multigrid, we need limit_level_difference_at_vertices
         // to make sure the transfer operators only need to consider two levels.
-      dealii::parallel::TriangulationBase<dim, spacedim>(
+      dealii::parallel::DistributedTriangulationBase<dim, spacedim>(
         mpi_communicator,
         (settings_ & construct_multigrid_hierarchy) ?
           static_cast<
@@ -2597,6 +2597,16 @@ namespace parallel
 
 
 
+    template <int dim, int spacedim>
+    bool
+    Triangulation<dim, spacedim>::is_multilevel_hierarchy_constructed() const
+    {
+      return settings &
+             Triangulation<dim, spacedim>::construct_multigrid_hierarchy;
+    }
+
+
+
     template <int dim, int spacedim>
     bool
     Triangulation<dim, spacedim>::has_hanging_nodes() const
@@ -4967,9 +4977,10 @@ namespace parallel
       const typename dealii::Triangulation<1, spacedim>::MeshSmoothing
         smooth_grid,
       const Settings /*settings*/)
-      : dealii::parallel::TriangulationBase<1, spacedim>(mpi_communicator,
-                                                         smooth_grid,
-                                                         false)
+      : dealii::parallel::DistributedTriangulationBase<1, spacedim>(
+          mpi_communicator,
+          smooth_grid,
+          false)
     {
       Assert(false, ExcNotImplemented());
     }
@@ -5104,6 +5115,16 @@ namespace parallel
       Assert(false, ExcNotImplemented());
     }
 
+
+
+    template <int spacedim>
+    bool
+    Triangulation<1, spacedim>::is_multilevel_hierarchy_constructed() const
+    {
+      Assert(false, ExcNotImplemented());
+      return false;
+    }
+
   } // namespace distributed
 } // namespace parallel
 
index 1fd06d882056381b823c7a46f71d311303e764fc..44252fa12f287f8aefa03f70cf322b3525c84e37 100644 (file)
@@ -392,6 +392,18 @@ namespace parallel
     return result;
   }
 
+  template <int dim, int spacedim>
+  DistributedTriangulationBase<dim, spacedim>::DistributedTriangulationBase(
+    MPI_Comm mpi_communicator,
+    const typename dealii::Triangulation<dim, spacedim>::MeshSmoothing
+               smooth_grid,
+    const bool check_for_distorted_cells)
+    : dealii::parallel::TriangulationBase<dim, spacedim>(
+        mpi_communicator,
+        smooth_grid,
+        check_for_distorted_cells)
+  {}
+
 } // end namespace parallel
 
 
index 3a6096f7ecbc166d2a2d7a04995271b9ff31ac37..861531c3928dc9229f2c285db0a3af9f37f73c2f 100644 (file)
@@ -28,5 +28,15 @@ for (deal_II_dimension : DIMENSIONS)
       template class TriangulationBase<deal_II_dimension,
                                        deal_II_dimension + 2>;
 #endif
+
+      template class DistributedTriangulationBase<deal_II_dimension>;
+#if deal_II_dimension < 3
+      template class DistributedTriangulationBase<deal_II_dimension,
+                                                  deal_II_dimension + 1>;
+#endif
+#if deal_II_dimension < 2
+      template class DistributedTriangulationBase<deal_II_dimension,
+                                                  deal_II_dimension + 2>;
+#endif
     \}
   }
index 4cb6d6228afe618155ab0cf29a4fa4428968720f..bc0c94a411bae7cce259ca6cec2f71f27fd449ff 100644 (file)
@@ -1446,7 +1446,7 @@ namespace internal
           std::vector<bool> include_vertex(
             dof_handler.get_triangulation().n_vertices(), false);
           if (dynamic_cast<
-                const parallel::distributed::Triangulation<dim, spacedim> *>(
+                const parallel::DistributedTriangulationBase<dim, spacedim> *>(
                 &dof_handler.get_triangulation()) != nullptr)
             for (const auto &cell : dof_handler.active_cell_iterators())
               if (cell->is_ghost())
@@ -3616,7 +3616,7 @@ namespace internal
         /**
          * level subdomain association. Similar to the above function only
          * for level meshes. This function assigns boundary dofs in
-         * the same way as parallel::distributed::Triangulation (proc with
+         * the same way as parallel::DistributedTriangulationBase (proc with
          * smallest index) instead of the coin flip method above.
          */
         template <class DoFHandlerType>
@@ -4366,7 +4366,7 @@ namespace internal
         template <int dim, int spacedim>
         void
         get_mg_dofindices_recursively(
-          const parallel::distributed::Triangulation<dim, spacedim> &tria,
+          const parallel::DistributedTriangulationBase<dim, spacedim> &tria,
           const typename dealii::internal::p4est::types<dim>::quadrant
             &p4est_cell,
           const typename DoFHandler<dim, spacedim>::level_cell_iterator
@@ -4420,7 +4420,7 @@ namespace internal
         template <int dim, int spacedim>
         void
         find_marked_mg_ghost_cells_recursively(
-          const typename parallel::distributed::Triangulation<dim, spacedim>
+          const typename parallel::DistributedTriangulationBase<dim, spacedim>
             &                tria,
           const unsigned int tree_index,
           const typename DoFHandler<dim, spacedim>::level_cell_iterator
@@ -4465,7 +4465,7 @@ namespace internal
         template <int dim, int spacedim>
         void
         set_mg_dofindices_recursively(
-          const parallel::distributed::Triangulation<dim, spacedim> &tria,
+          const parallel::DistributedTriangulationBase<dim, spacedim> &tria,
           const typename dealii::internal::p4est::types<dim>::quadrant
             &p4est_cell,
           const typename DoFHandler<dim, spacedim>::level_cell_iterator
@@ -4536,7 +4536,7 @@ namespace internal
         template <int dim, int spacedim, class DoFHandlerType>
         void
         communicate_mg_ghost_cells(
-          const typename parallel::distributed::Triangulation<dim, spacedim>
+          const typename parallel::DistributedTriangulationBase<dim, spacedim>
             &             tria,
           DoFHandlerType &dof_handler)
         {
@@ -4961,7 +4961,7 @@ namespace internal
           // barrier is negligible compared to everything else we do
           // here
           if (const auto *triangulation = dynamic_cast<
-                const parallel::distributed::Triangulation<dim, spacedim> *>(
+                const parallel::DistributedTriangulationBase<dim, spacedim> *>(
                 &dof_handler.get_triangulation()))
             {
               const int ierr = MPI_Barrier(triangulation->get_communicator());
@@ -5004,10 +5004,10 @@ namespace internal
         const unsigned int dim      = DoFHandlerType::dimension;
         const unsigned int spacedim = DoFHandlerType::space_dimension;
 
-        parallel::distributed::Triangulation<dim, spacedim> *triangulation =
-          (dynamic_cast<parallel::distributed::Triangulation<dim, spacedim> *>(
-            const_cast<dealii::Triangulation<dim, spacedim> *>(
-              &dof_handler->get_triangulation())));
+        parallel::DistributedTriangulationBase<dim, spacedim> *triangulation =
+          (dynamic_cast<parallel::DistributedTriangulationBase<dim, spacedim>
+                          *>(const_cast<dealii::Triangulation<dim, spacedim> *>(
+            &dof_handler->get_triangulation())));
         Assert(triangulation != nullptr, ExcInternalError());
 
         const types::subdomain_id subdomain_id =
@@ -5241,15 +5241,13 @@ namespace internal
         const unsigned int dim      = DoFHandlerType::dimension;
         const unsigned int spacedim = DoFHandlerType::space_dimension;
 
-        parallel::distributed::Triangulation<dim, spacedim> *triangulation =
-          (dynamic_cast<parallel::distributed::Triangulation<dim, spacedim> *>(
-            const_cast<dealii::Triangulation<dim, spacedim> *>(
-              &dof_handler->get_triangulation())));
+        parallel::DistributedTriangulationBase<dim, spacedim> *triangulation =
+          (dynamic_cast<parallel::DistributedTriangulationBase<dim, spacedim>
+                          *>(const_cast<dealii::Triangulation<dim, spacedim> *>(
+            &dof_handler->get_triangulation())));
         Assert(triangulation != nullptr, ExcInternalError());
 
-        AssertThrow((triangulation->settings &
-                     parallel::distributed::Triangulation<dim, spacedim>::
-                       construct_multigrid_hierarchy),
+        AssertThrow((triangulation->is_multilevel_hierarchy_constructed()),
                     ExcMessage(
                       "Multigrid DoFs can only be distributed on a parallel "
                       "Triangulation if the flag construct_multigrid_hierarchy "
@@ -5470,10 +5468,10 @@ namespace internal
         const unsigned int dim      = DoFHandlerType::dimension;
         const unsigned int spacedim = DoFHandlerType::space_dimension;
 
-        parallel::distributed::Triangulation<dim, spacedim> *triangulation =
-          (dynamic_cast<parallel::distributed::Triangulation<dim, spacedim> *>(
-            const_cast<dealii::Triangulation<dim, spacedim> *>(
-              &dof_handler->get_triangulation())));
+        parallel::DistributedTriangulationBase<dim, spacedim> *triangulation =
+          (dynamic_cast<parallel::DistributedTriangulationBase<dim, spacedim>
+                          *>(const_cast<dealii::Triangulation<dim, spacedim> *>(
+            &dof_handler->get_triangulation())));
         Assert(triangulation != nullptr, ExcInternalError());
 
 

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.