]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Exchange ghost uses vertices w ghost_neighbors
authorBruno Blais <blais.bruno@gmail.com>
Mon, 21 Nov 2022 01:41:51 +0000 (20:41 -0500)
committerBruno Blais <blais.bruno@gmail.com>
Mon, 21 Nov 2022 01:41:51 +0000 (20:41 -0500)
include/deal.II/grid/grid_tools_cache.h
include/deal.II/grid/grid_tools_cache_update_flags.h
source/grid/grid_tools_cache.cc
source/particles/particle_handler.cc

index 6d0a3e4392f4e2072221437c0c4f1f180ba3b586..7b44dbda1df0626bd7885a2ecd85067b0f7a22bf 100644 (file)
@@ -75,7 +75,7 @@ namespace GridTools
      * @param mapping The mapping to use when computing cached objects
      */
     Cache(const Triangulation<dim, spacedim> &tria,
-          const Mapping<dim, spacedim> &      mapping =
+          const Mapping<dim, spacedim>       &mapping =
             (ReferenceCells::get_hypercube<dim>()
 #ifndef _MSC_VER
                .template get_default_linear_mapping<dim, spacedim>()
@@ -169,22 +169,11 @@ namespace GridTools
     get_vertex_to_neighbor_subdomain() const;
 
     /**
-     * Returns the std::map of std::vector of unsigned integer containing
-     * coinciding vertices labeled by an arbitrary element from them. This
-     * feature is used to identify cells which are connected by periodic
-     * boundaries.
+     * Return a map that, for each vertex, lists all the processes whose
+     * subdomains are adjacent to that vertex.
      */
-    const std::map<unsigned int, std::vector<unsigned int>> &
-    get_coinciding_vertex_groups() const;
-
-    /**
-     * Returns the std::map of unsigned integer containing
-     * the label of a group of coinciding vertices. This
-     * feature is used to identify cells which are connected by periodic
-     * boundaries.
-     */
-    const std::map<unsigned int, unsigned int> &
-    get_vertex_to_coinciding_vertex_group() const;
+    const std::map<unsigned int, std::set<types::subdomain_id>> &
+    get_vertices_with_ghost_neighbors() const;
 
     /**
      * Return a reference to the stored triangulation.
@@ -310,21 +299,12 @@ namespace GridTools
      */
     mutable std::vector<std::set<unsigned int>> vertex_to_neighbor_subdomain;
 
-    /**
-     * Store an std::map of std::vector of unsigned integer containing
-     * coinciding vertices labeled by an arbitrary element from them.
-     */
-    mutable std::map<unsigned int, std::vector<unsigned int>>
-      coinciding_vertex_groups;
-
     /**
      * Store an std::map of unsigned integer containing
-     * the label of a group of coinciding vertices.
-     * This std::map is generally used in combination
-     * with the above std::map coinciding_vertex_groups.
+     * the set of subdomains connected to each vertices
      */
-    mutable std::map<unsigned int, unsigned int>
-      vertex_to_coinciding_vertex_group;
+    mutable std::map<unsigned int, std::set<dealii::types::subdomain_id>>
+      vertices_with_ghost_neighbors;
 
     /**
      * Storage for the status of the triangulation signal.
index 381c3543b1e30dce5ea04698cc8083ff210db921..9e5750c56a18a6fdba158adb90622eeb530ca96b 100644 (file)
@@ -88,7 +88,7 @@ namespace GridTools
      * Update the coinciding vertex groups as well as the
      * vertex_to_coinciding_vertex_group
      */
-    update_collection_of_coinciding_vertices = 0x200,
+    update_vertex_with_ghost_neighbors = 0x200,
 
     /**
      * Update all objects.
index 0b8762b4a859f6ace9edd688db0ac4e8d0fd3977..f808f2889c1100c8557194d69dc3d553d91d2a65 100644 (file)
@@ -26,7 +26,7 @@ namespace GridTools
 {
   template <int dim, int spacedim>
   Cache<dim, spacedim>::Cache(const Triangulation<dim, spacedim> &tria,
-                              const Mapping<dim, spacedim> &      mapping)
+                              const Mapping<dim, spacedim>       &mapping)
     : update_flags(update_all)
     , tria(&tria)
     , mapping(&mapping)
@@ -219,33 +219,18 @@ namespace GridTools
   }
 
   template <int dim, int spacedim>
-  const std::map<unsigned int, std::vector<unsigned int>> &
-  Cache<dim, spacedim>::get_coinciding_vertex_groups() const
+  const std::map<unsigned int, std::set<types::subdomain_id>> &
+  Cache<dim, spacedim>::get_vertices_with_ghost_neighbors() const
   {
-    if (update_flags & update_collection_of_coinciding_vertices)
+    if (update_flags & update_vertex_with_ghost_neighbors)
       {
-        GridTools::collect_coinciding_vertices(
-          *tria, coinciding_vertex_groups, vertex_to_coinciding_vertex_group);
+        vertices_with_ghost_neighbors =
+          GridTools::compute_vertices_with_ghost_neighbors(*tria);
 
-        update_flags = update_flags & ~update_collection_of_coinciding_vertices;
+        update_flags = update_flags & ~update_vertex_with_ghost_neighbors;
       }
 
-    return coinciding_vertex_groups;
-  }
-
-  template <int dim, int spacedim>
-  const std::map<unsigned int, unsigned int> &
-  Cache<dim, spacedim>::get_vertex_to_coinciding_vertex_group() const
-  {
-    if (update_flags & update_collection_of_coinciding_vertices)
-      {
-        GridTools::collect_coinciding_vertices(
-          *tria, coinciding_vertex_groups, vertex_to_coinciding_vertex_group);
-
-        update_flags = update_flags & ~update_collection_of_coinciding_vertices;
-      }
-
-    return vertex_to_coinciding_vertex_group;
+    return vertices_with_ghost_neighbors;
   }
 
 #include "grid_tools_cache.inst"
index a1c710ccf82cb508c41c1a30aff5752257abed72..320dda647fc74bbcb6f094367801ea2ac505669b 100644 (file)
@@ -77,7 +77,7 @@ namespace Particles
   template <int dim, int spacedim>
   ParticleHandler<dim, spacedim>::ParticleHandler(
     const Triangulation<dim, spacedim> &triangulation,
-    const Mapping<dim, spacedim> &      mapping,
+    const Mapping<dim, spacedim>       &mapping,
     const unsigned int                  n_properties)
     : triangulation(&triangulation, typeid(*this).name())
     , mapping(&mapping, typeid(*this).name())
@@ -117,7 +117,7 @@ namespace Particles
   void
   ParticleHandler<dim, spacedim>::initialize(
     const Triangulation<dim, spacedim> &new_triangulation,
-    const Mapping<dim, spacedim> &      new_mapping,
+    const Mapping<dim, spacedim>       &new_mapping,
     const unsigned int                  n_properties)
   {
     clear();
@@ -577,7 +577,7 @@ namespace Particles
   template <int dim, int spacedim>
   typename ParticleHandler<dim, spacedim>::particle_iterator
   ParticleHandler<dim, spacedim>::insert_particle(
-    const Particle<dim, spacedim> &                                    particle,
+    const Particle<dim, spacedim>                                     &particle,
     const typename Triangulation<dim, spacedim>::active_cell_iterator &cell)
   {
     return insert_particle(particle.get_location(),
@@ -623,7 +623,7 @@ namespace Particles
   template <int dim, int spacedim>
   typename ParticleHandler<dim, spacedim>::particle_iterator
   ParticleHandler<dim, spacedim>::insert_particle(
-    const void *&                                                      data,
+    const void                                                       *&data,
     const typename Triangulation<dim, spacedim>::active_cell_iterator &cell)
   {
     Assert(triangulation != nullptr, ExcInternalError());
@@ -648,8 +648,8 @@ namespace Particles
   template <int dim, int spacedim>
   typename ParticleHandler<dim, spacedim>::particle_iterator
   ParticleHandler<dim, spacedim>::insert_particle(
-    const Point<spacedim> &     position,
-    const Point<dim> &          reference_position,
+    const Point<spacedim>      &position,
+    const Point<dim>           &reference_position,
     const types::particle_index particle_index,
     const typename Triangulation<dim, spacedim>::active_cell_iterator &cell,
     const ArrayView<const double> &properties)
@@ -764,8 +764,8 @@ namespace Particles
   ParticleHandler<dim, spacedim>::insert_global_particles(
     const std::vector<Point<spacedim>> &positions,
     const std::vector<std::vector<BoundingBox<spacedim>>>
-      &                                       global_bounding_boxes,
-    const std::vector<std::vector<double>> &  properties,
+                                             &global_bounding_boxes,
+    const std::vector<std::vector<double>>   &properties,
     const std::vector<types::particle_index> &ids)
   {
     if (!properties.empty())
@@ -1190,7 +1190,7 @@ namespace Particles
     compare_particle_association(
       const unsigned int                 a,
       const unsigned int                 b,
-      const Tensor<1, dim> &             particle_direction,
+      const Tensor<1, dim>              &particle_direction,
       const std::vector<Tensor<1, dim>> &center_directions)
     {
       const double scalar_product_a = center_directions[a] * particle_direction;
@@ -1545,14 +1545,11 @@ namespace Particles
     // In the case of a parallel simulation with periodic boundary conditions
     // the vertices associated with periodic boundaries are not directly
     // connected to the ghost cells but they are connected to the ghost cells
-    // through their coinciding vertices. We gather the information about
-    // the coinciding vertices through the grid cache.
-    const std::map<unsigned int, std::vector<unsigned int>>
-      &coinciding_vertex_groups =
-        triangulation_cache->get_coinciding_vertex_groups();
-    const std::map<unsigned int, unsigned int>
-      &vertex_to_coinciding_vertex_group =
-        triangulation_cache->get_vertex_to_coinciding_vertex_group();
+    // through their coinciding vertices. We gather this information using the
+    // vertices_with_ghost_neighbors map
+    const std::map<unsigned int, std::set<types::subdomain_id>>
+      &vertices_with_ghost_neighbors =
+        triangulation_cache->get_vertices_with_ghost_neighbors();
 
     const std::set<types::subdomain_id> ghost_owners =
       parallel_triangulation->ghost_owners();
@@ -1570,29 +1567,14 @@ namespace Particles
             std::set<unsigned int> cell_to_neighbor_subdomain;
             for (const unsigned int v : cell->vertex_indices())
               {
-                cell_to_neighbor_subdomain.insert(
-                  vertex_to_neighbor_subdomain[cell->vertex_index(v)].begin(),
-                  vertex_to_neighbor_subdomain[cell->vertex_index(v)].end());
-
-                // If the vertex has one or more coinciding vertices, then
-                // all of these vertices, including the current vertex being
-                // looped over are contained inside the
-                // vertex_to_coinciding_vertex_group. Consequently, we loop over
-                // the full content of the coinciding vertex group which will
-                // include the current vertex. When we reach the current vertex,
-                // we skip it to avoid repetition.
-                if (vertex_to_coinciding_vertex_group.find(cell->vertex_index(
-                      v)) != vertex_to_coinciding_vertex_group.end())
+                if (vertices_with_ghost_neighbors.find(cell->vertex_index(v)) !=
+                    vertices_with_ghost_neighbors.end())
                   {
-                    for (const unsigned int cv : coinciding_vertex_groups.at(
-                           vertex_to_coinciding_vertex_group.at(
-                             cell->vertex_index(v))))
-                      {
-                        if (cv != cell->vertex_index(v))
-                          cell_to_neighbor_subdomain.insert(
-                            vertex_to_neighbor_subdomain[cv].begin(),
-                            vertex_to_neighbor_subdomain[cv].end());
-                      }
+                    cell_to_neighbor_subdomain.insert(
+                      vertices_with_ghost_neighbors.at(cell->vertex_index(v))
+                        .begin(),
+                      vertices_with_ghost_neighbors.at(cell->vertex_index(v))
+                        .end());
                   }
               }
 
@@ -1669,7 +1651,7 @@ namespace Particles
     const std::map<
       types::subdomain_id,
       std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>>
-      &        send_cells,
+              &send_cells,
     const bool build_cache)
   {
     Assert(triangulation != nullptr, ExcInternalError());
@@ -2073,7 +2055,7 @@ namespace Particles
   template <int dim, int spacedim>
   void
   ParticleHandler<dim, spacedim>::register_additional_store_load_functions(
-    const std::function<std::size_t()> &                            size_callb,
+    const std::function<std::size_t()>                             &size_callb,
     const std::function<void *(const particle_iterator &, void *)> &store_callb,
     const std::function<const void *(const particle_iterator &, const void *)>
       &load_callb)
@@ -2203,7 +2185,7 @@ namespace Particles
     const auto callback_function =
       [this](
         const typename Triangulation<dim, spacedim>::cell_iterator
-          &                                                     cell_iterator,
+                                                               &cell_iterator,
         const typename Triangulation<dim, spacedim>::CellStatus cell_status) {
         return this->pack_callback(cell_iterator, cell_status);
       };
@@ -2360,7 +2342,7 @@ namespace Particles
   template <int dim, int spacedim>
   void
   ParticleHandler<dim, spacedim>::unpack_callback(
-    const typename Triangulation<dim, spacedim>::cell_iterator &    cell,
+    const typename Triangulation<dim, spacedim>::cell_iterator     &cell,
     const typename Triangulation<dim, spacedim>::CellStatus         status,
     const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
   {
@@ -2394,14 +2376,12 @@ namespace Particles
     // now update particle storage location and properties if necessary
     switch (status)
       {
-        case parallel::TriangulationBase<dim, spacedim>::CELL_PERSIST:
-          {
+          case parallel::TriangulationBase<dim, spacedim>::CELL_PERSIST: {
             // all particles are correctly inserted
           }
           break;
 
-        case parallel::TriangulationBase<dim, spacedim>::CELL_COARSEN:
-          {
+          case parallel::TriangulationBase<dim, spacedim>::CELL_COARSEN: {
             // all particles are in correct cell, but their reference location
             // has changed
             for (auto &particle : loaded_particles_on_cell)
@@ -2414,8 +2394,7 @@ namespace Particles
           }
           break;
 
-        case parallel::TriangulationBase<dim, spacedim>::CELL_REFINE:
-          {
+          case parallel::TriangulationBase<dim, spacedim>::CELL_REFINE: {
             // we need to find the correct child to store the particles and
             // their reference location has changed
             typename particle_container::iterator &cache =

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.