]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Triangulation: implement global_*_cell_index_partitioner() 14270/head
authorPeter Munch <peterrmuench@gmail.com>
Wed, 14 Sep 2022 07:08:27 +0000 (09:08 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Wed, 14 Sep 2022 19:41:33 +0000 (21:41 +0200)
include/deal.II/distributed/tria_base.h
include/deal.II/grid/tria.h
source/grid/tria.cc

index 011c2d32298ec44333a95204f466ab01714161cc..012de34920c8a5bad55a754d2a43a8419efd175c 100644 (file)
@@ -202,21 +202,12 @@ namespace parallel
     const std::set<types::subdomain_id> &
     level_ghost_owners() const;
 
-    /**
-     * Return partitioner for the global indices of the cells on the active
-     * level of the triangulation, which is returned by the function
-     * CellAccessor::global_active_cell_index().
-     */
     const std::weak_ptr<const Utilities::MPI::Partitioner>
-    global_active_cell_index_partitioner() const;
+    global_active_cell_index_partitioner() const override;
 
-    /**
-     * Return partitioner for the global indices of the cells on the given @p
-     * level of the triangulation, which is returned by the function
-     * CellAccessor::global_level_cell_index().
-     */
     const std::weak_ptr<const Utilities::MPI::Partitioner>
-    global_level_cell_index_partitioner(const unsigned int level) const;
+    global_level_cell_index_partitioner(
+      const unsigned int level) const override;
 
     /**
      * @copydoc dealii::Triangulation::get_boundary_ids()
index ae36c0c9b86bcc7317a0f99b4d9ec2d61ae19f22..3f792c8e6aadf04469e1cfe38449c5d4e8abf770 100644 (file)
@@ -21,6 +21,7 @@
 
 #include <deal.II/base/geometry_info.h>
 #include <deal.II/base/iterator_range.h>
+#include <deal.II/base/partitioner.h>
 #include <deal.II/base/point.h>
 #include <deal.II/base/smartpointer.h>
 #include <deal.II/base/subscriptor.h>
@@ -176,6 +177,18 @@ namespace internal
        */
       std::vector<unsigned int> n_active_lines_level;
 
+      /**
+       * Partitioner for the global active cell indices.
+       */
+      std::shared_ptr<const Utilities::MPI::Partitioner>
+        active_cell_index_partitioner;
+
+      /**
+       * Partitioner for the global level cell indices for each level.
+       */
+      std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
+        level_cell_index_partitioners;
+
       /**
        * Constructor. Set values to zero by default.
        */
@@ -1615,6 +1628,22 @@ public:
   virtual MPI_Comm
   get_communicator() const;
 
+  /**
+   * Return the partitioner for the global indices of the cells on the active
+   * level of the triangulation, which is returned by the function
+   * CellAccessor::global_active_cell_index().
+   */
+  virtual const std::weak_ptr<const Utilities::MPI::Partitioner>
+  global_active_cell_index_partitioner() const;
+
+  /**
+   * Return the partitioner for the global indices of the cells on the given @p
+   * level of the triangulation, which is returned by the function
+   * CellAccessor::global_level_cell_index().
+   */
+  virtual const std::weak_ptr<const Utilities::MPI::Partitioner>
+  global_level_cell_index_partitioner(const unsigned int level) const;
+
   /**
    * Set the mesh smoothing to @p mesh_smoothing. This overrides the
    * MeshSmoothing given to the constructor. It is allowed to call this
index 44abfd59e7af6a7d9317e364911e686fe2bc1343..87f10c5d230a78445cc3492a5f629f3d6c42fea4 100644 (file)
@@ -1757,7 +1757,7 @@ namespace internal
        */
       template <int dim, int spacedim>
       static void
-      compute_number_cache(
+      compute_number_cache_dim(
         const Triangulation<dim, spacedim> &                   triangulation,
         const unsigned int                                     level_objects,
         internal::TriangulationImplementation::NumberCache<1> &number_cache)
@@ -1845,7 +1845,7 @@ namespace internal
        */
       template <int dim, int spacedim>
       static void
-      compute_number_cache(
+      compute_number_cache_dim(
         const Triangulation<dim, spacedim> &                   triangulation,
         const unsigned int                                     level_objects,
         internal::TriangulationImplementation::NumberCache<2> &number_cache)
@@ -1858,7 +1858,7 @@ namespace internal
             void (*)(const Triangulation<dim, spacedim> &,
                      const unsigned int,
                      internal::TriangulationImplementation::NumberCache<1> &)>(
-            &compute_number_cache<dim, spacedim>),
+            &compute_number_cache_dim<dim, spacedim>),
           triangulation,
           level_objects,
           static_cast<internal::TriangulationImplementation::NumberCache<1> &>(
@@ -1952,7 +1952,7 @@ namespace internal
        */
       template <int dim, int spacedim>
       static void
-      compute_number_cache(
+      compute_number_cache_dim(
         const Triangulation<dim, spacedim> &                   triangulation,
         const unsigned int                                     level_objects,
         internal::TriangulationImplementation::NumberCache<3> &number_cache)
@@ -1965,7 +1965,7 @@ namespace internal
             void (*)(const Triangulation<dim, spacedim> &,
                      const unsigned int,
                      internal::TriangulationImplementation::NumberCache<2> &)>(
-            &compute_number_cache<dim, spacedim>),
+            &compute_number_cache_dim<dim, spacedim>),
           triangulation,
           level_objects,
           static_cast<internal::TriangulationImplementation::NumberCache<2> &>(
@@ -2043,6 +2043,27 @@ namespace internal
       }
 
 
+      template <int dim, int spacedim>
+      static void
+      compute_number_cache(
+        const Triangulation<dim, spacedim> &                     triangulation,
+        const unsigned int                                       level_objects,
+        internal::TriangulationImplementation::NumberCache<dim> &number_cache)
+      {
+        compute_number_cache_dim(triangulation, level_objects, number_cache);
+
+        number_cache.active_cell_index_partitioner =
+          std::make_shared<const Utilities::MPI::Partitioner>(
+            triangulation.n_active_cells());
+
+        number_cache.level_cell_index_partitioners.resize(
+          triangulation.n_levels());
+        for (unsigned int level = 0; level < triangulation.n_levels(); ++level)
+          number_cache.level_cell_index_partitioners[level] =
+            std::make_shared<const Utilities::MPI::Partitioner>(
+              triangulation.n_cells(level));
+      }
+
 
       template <int spacedim>
       static void
@@ -11207,6 +11228,27 @@ Triangulation<dim, spacedim>::get_communicator() const
 
 
 
+template <int dim, int spacedim>
+const std::weak_ptr<const Utilities::MPI::Partitioner>
+Triangulation<dim, spacedim>::global_active_cell_index_partitioner() const
+{
+  return number_cache.active_cell_index_partitioner;
+}
+
+
+
+template <int dim, int spacedim>
+const std::weak_ptr<const Utilities::MPI::Partitioner>
+Triangulation<dim, spacedim>::global_level_cell_index_partitioner(
+  const unsigned int level) const
+{
+  AssertIndexRange(level, this->n_levels());
+
+  return number_cache.level_cell_index_partitioners[level];
+}
+
+
+
 template <int dim, int spacedim>
 void
 Triangulation<dim, spacedim>::set_mesh_smoothing(

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.