]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add Triangulation::n_global_coarse_cells() 12502/head
authorPeter Munch <peterrmuench@gmail.com>
Sat, 26 Jun 2021 15:30:16 +0000 (17:30 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 28 Jun 2021 17:40:50 +0000 (19:40 +0200)
doc/news/changes/minor/20210626Munch [new file with mode: 0644]
include/deal.II/distributed/fully_distributed_tria.h
include/deal.II/distributed/tria_base.h
include/deal.II/grid/tria.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
source/distributed/fully_distributed_tria.cc
source/distributed/repartitioning_policy_tools.cc
source/distributed/tria_base.cc
source/grid/tria.cc

diff --git a/doc/news/changes/minor/20210626Munch b/doc/news/changes/minor/20210626Munch
new file mode 100644 (file)
index 0000000..c57fdc5
--- /dev/null
@@ -0,0 +1,6 @@
+New: The new function Triangulation::n_global_coarse_cells() allows to
+query the global number of coarse cells. This function is particularly useful
+for parallel::fullydistributed::Triangulation, since in this case the coarse 
+cells might differ on each process.
+<br>
+(Peter Munch, 2021/06/26)
index 27414c756ef033c085c7e11b1778d7feb79c43b2..326dec1120c66e2d0344591b8050fe86c38cb8c6 100644 (file)
@@ -278,6 +278,9 @@ namespace parallel
       virtual void
       update_cell_relations() override;
 
+      virtual void
+      update_number_cache() override;
+
       /**
        * store the Settings.
        */
index 6e34a49e7a84e91c0659b99142337abcb75a536a..076a4e6e06daf24c1713df72896d0a8fc2f3ef6f 100644 (file)
@@ -289,6 +289,9 @@ namespace parallel
     communicate_locally_moved_vertices(
       const std::vector<bool> &vertex_locally_moved);
 
+    virtual types::coarse_cell_id
+    n_global_coarse_cells() const override;
+
   protected:
     /**
      * MPI communicator to be used for the triangulation. We create a unique
@@ -318,22 +321,31 @@ namespace parallel
        * Number of locally owned active cells of this MPI rank.
        */
       unsigned int n_locally_owned_active_cells;
+
       /**
        * The total number of active cells (sum of @p
        * n_locally_owned_active_cells).
        */
       types::global_cell_index n_global_active_cells;
+
+      /**
+       * Number of global coarse cells.
+       */
+      types::coarse_cell_id number_of_global_coarse_cells;
+
       /**
        * The global number of levels computed as the maximum number of levels
        * taken over all MPI ranks, so <tt>n_levels()<=n_global_levels =
        * max(n_levels() on proc i)</tt>.
        */
       unsigned int n_global_levels;
+
       /**
        * A set containing the subdomain_id (MPI rank) of the owners of the
        * ghost cells on this processor.
        */
       std::set<types::subdomain_id> ghost_owners;
+
       /**
        * A set containing the MPI ranks of the owners of the level ghost cells
        * on this processor (for all levels).
index b954f00508ce6e227e036f23f239e23958291588..265b70d0ae650be17ffaa2b90b75d64f74327c95 100644 (file)
@@ -3082,6 +3082,13 @@ public:
   unsigned int
   n_active_cells(const unsigned int level) const;
 
+  /**
+   * Return the total number of coarse cells. If the coarse mesh is replicated
+   * on each process, this simply returns <tt>n_cells(0)</tt>.
+   */
+  virtual types::coarse_cell_id
+  n_global_coarse_cells() const;
+
   /**
    * Return the total number of used faces, active or not.  In 2D, the result
    * equals n_lines(), in 3D it equals n_quads(), while in 1D it equals
index 58ef8e919215fff183e8bf51223e0f1f245a7a32..cbb38d2de2cbf91bf7585a27140ac34d377550b3 100644 (file)
@@ -551,12 +551,14 @@ namespace internal
       , mg_level_fine(mg_level_fine)
       , communicator(
           mesh_fine.get_communicator() /*TODO: fix for different comms*/)
-      , cell_id_translator(n_coarse_cells(mesh_fine),
-                           n_global_levels(mesh_fine))
+      , cell_id_translator(
+          mesh_fine.get_triangulation().n_global_coarse_cells(),
+          mesh_fine.get_triangulation().n_global_levels())
     {
-      AssertDimension(n_coarse_cells(mesh_fine), n_coarse_cells(mesh_coarse));
-      AssertIndexRange(n_global_levels(mesh_coarse),
-                       n_global_levels(mesh_fine) + 1);
+      AssertDimension(mesh_fine.get_triangulation().n_global_coarse_cells(),
+                      mesh_coarse.get_triangulation().n_global_coarse_cells());
+      AssertIndexRange(mesh_coarse.get_triangulation().n_global_levels(),
+                       mesh_fine.get_triangulation().n_global_levels() + 1);
     }
 
     void
@@ -969,25 +971,6 @@ namespace internal
     std::map<unsigned int,
              std::pair<unsigned int, std::vector<types::global_dof_index>>>
       map;
-
-    static unsigned int
-    n_coarse_cells(const DoFHandler<dim> &mesh)
-    {
-      types::coarse_cell_id n_coarse_cells = 0;
-
-      for (const auto &cell : mesh.get_triangulation().active_cell_iterators())
-        if (!cell->is_artificial())
-          n_coarse_cells =
-            std::max(n_coarse_cells, cell->id().get_coarse_cell_id());
-
-      return Utilities::MPI::max(n_coarse_cells, mesh.get_communicator()) + 1;
-    }
-
-    static unsigned int
-    n_global_levels(const DoFHandler<dim> &mesh)
-    {
-      return mesh.get_triangulation().n_global_levels();
-    }
   };
 
   template <int dim>
index 95ed52fcceb948a8a6b1fe582cdd4387471fd5ad..22033a0afaf5cae79594e4fb004d81292e142eac 100644 (file)
@@ -701,6 +701,32 @@ namespace parallel
     }
 
 
+
+    template <int dim, int spacedim>
+    void
+    Triangulation<dim, spacedim>::update_number_cache()
+    {
+      dealii::parallel::TriangulationBase<dim, spacedim>::update_number_cache();
+
+      // additionally update the number of global coarse cells
+      types::coarse_cell_id number_of_global_coarse_cells = 0;
+
+      for (const auto &cell : this->active_cell_iterators())
+        if (!cell->is_artificial())
+          number_of_global_coarse_cells =
+            std::max(number_of_global_coarse_cells,
+                     cell->id().get_coarse_cell_id());
+
+      number_of_global_coarse_cells =
+        Utilities::MPI::max(number_of_global_coarse_cells,
+                            this->mpi_communicator) +
+        1;
+
+      this->number_cache.number_of_global_coarse_cells =
+        number_of_global_coarse_cells;
+    }
+
+
   } // namespace fullydistributed
 } // namespace parallel
 
index a7e35036f42772e8ff6a3a16848b3faec8c7c6af..31e4739d156ab5aa23235b5917f06f58253f85ad 100644 (file)
@@ -28,31 +28,6 @@ namespace RepartitioningPolicyTools
 {
   namespace
   {
-    template <typename MeshType>
-    unsigned int
-    compute_n_coarse_cells(const MeshType &mesh)
-    {
-      types::coarse_cell_id n_coarse_cells = 0;
-
-      for (const auto &cell :
-           mesh.get_triangulation().cell_iterators_on_level(0))
-        n_coarse_cells =
-          std::max(n_coarse_cells, cell->id().get_coarse_cell_id());
-
-      return Utilities::MPI::max(n_coarse_cells, mesh.get_communicator()) + 1;
-    }
-
-
-
-    template <typename MeshType>
-    unsigned int
-    compute_n_global_levels(const MeshType &mesh)
-    {
-      return mesh.get_triangulation().n_global_levels();
-    }
-
-
-
     template <int dim, int spacedim>
     void
     add_indices_recursevly_for_first_child_policy(
@@ -83,8 +58,8 @@ namespace RepartitioningPolicyTools
   template <int dim, int spacedim>
   FirstChildPolicy<dim, spacedim>::FirstChildPolicy(
     const Triangulation<dim, spacedim> &tria_fine)
-    : n_coarse_cells(compute_n_coarse_cells(tria_fine))
-    , n_global_levels(compute_n_global_levels(tria_fine))
+    : n_coarse_cells(tria_fine.n_global_coarse_cells())
+    , n_global_levels(tria_fine.n_global_levels())
   {
     Assert(
       tria_fine.all_reference_cells_are_hyper_cube(),
index 864e982435c4f32fbd5befe42ae3ef9eee21de3c..0035491eb92aef78cf519820d7dda7273e5d4334 100644 (file)
@@ -111,6 +111,7 @@ namespace parallel
   template <int dim, int spacedim>
   TriangulationBase<dim, spacedim>::NumberCache::NumberCache()
     : n_locally_owned_active_cells(0)
+    , number_of_global_coarse_cells(0)
     , n_global_levels(0)
   {}
 
@@ -273,6 +274,8 @@ namespace parallel
                ExcInternalError());
       }
 
+    this->number_cache.number_of_global_coarse_cells = this->n_cells(0);
+
     // reset global cell ids
     this->reset_global_cell_indices();
   }
@@ -649,6 +652,15 @@ namespace parallel
 
 
 
+  template <int dim, int spacedim>
+  types::coarse_cell_id
+  TriangulationBase<dim, spacedim>::n_global_coarse_cells() const
+  {
+    return number_cache.number_of_global_coarse_cells;
+  }
+
+
+
   template <int dim, int spacedim>
   DistributedTriangulationBase<dim, spacedim>::DistributedTriangulationBase(
     const MPI_Comm &mpi_communicator,
index b5c3724acad752ad07dac1b2c6e0535c1175136a..079b622b6c1bd5a4ab766ab7036d88841ce978c4 100644 (file)
@@ -13830,7 +13830,12 @@ Triangulation<dim, spacedim>::n_global_active_cells() const
   return n_active_cells();
 }
 
-
+template <int dim, int spacedim>
+types::coarse_cell_id
+Triangulation<dim, spacedim>::n_global_coarse_cells() const
+{
+  return n_cells(0);
+}
 
 template <int dim, int spacedim>
 unsigned 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.