]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added vertex_to_neighbhor_subdomain to cache and use it therein the
authorBruno <blais.bruno@gmail.com>
Sun, 25 Oct 2020 18:39:22 +0000 (14:39 -0400)
committerBruno Blais <blais.bruno@gmail.com>
Wed, 28 Oct 2020 12:31:15 +0000 (08:31 -0400)
particle_handler

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 26140d8097a4a8c2726157394b902b12698da143..f61b1f3403fe451bf0024ce7354bf32a4de1279c 100644 (file)
@@ -152,6 +152,15 @@ namespace GridTools
                 typename Triangulation<dim, spacedim>::active_cell_iterator>> &
     get_locally_owned_cell_bounding_boxes_rtree() const;
 
+
+    /** Returns the vector of set of integer containing the subdomain id
+     * to which each vertex is connected to. This feature is used extensively
+     * in the particle_handler to detect on which processors ghost particles
+     * must be built
+     */
+    const std::vector<std::set<unsigned int>> &
+    get_vertex_to_neighbor_subdomain() const;
+
     /**
      * Return a reference to the stored triangulation.
      */
@@ -269,6 +278,13 @@ namespace GridTools
                 typename Triangulation<dim, spacedim>::active_cell_iterator>>
       locally_owned_cell_bounding_boxes_rtree;
 
+
+    /**
+     * Store an std::vector of std::set of integer containing the id of all
+     * subdomain to which a vertex is connected to
+     */
+    mutable std::vector<std::set<unsigned int>> vertex_to_neighbor_subdomain;
+
     /**
      * Storage for the status of the triangulation signal.
      */
index 2a89e3236c52e9d53434622b78e804da2ad58fd5..c8d3e4c05b91fd97808fcebae23f98881e1bf235 100644 (file)
@@ -79,6 +79,11 @@ namespace GridTools
      */
     update_locally_owned_cell_bounding_boxes_rtree = 0x080,
 
+    /**
+     * Update vertex to neighbhor subdomain
+     */
+    update_vertex_to_neighbor_subdomain = 0x100,
+
     /**
      * Update all objects.
      */
@@ -163,8 +168,8 @@ namespace GridTools
    *
    * @ref CacheUpdateFlags
    */
-  inline CacheUpdateFlags operator&(const CacheUpdateFlags f1,
-                                    const CacheUpdateFlags f2)
+  inline CacheUpdateFlags
+  operator&(const CacheUpdateFlags f1, const CacheUpdateFlags f2)
   {
     return static_cast<CacheUpdateFlags>(static_cast<unsigned int>(f1) &
                                          static_cast<unsigned int>(f2));
index 9cd4a8704bc3eb3ee2049c4dee4b4e5a257108ac..46986e3ca40ad2b84a17f17329a078efa507d3a1 100644 (file)
@@ -199,6 +199,25 @@ namespace GridTools
     return covering_rtree[level];
   }
 
+  template <int dim, int spacedim>
+  const std::vector<std::set<unsigned int>> &
+  Cache<dim, spacedim>::get_vertex_to_neighbor_subdomain() const
+  {
+    if (update_flags & update_vertex_to_neighbor_subdomain)
+      {
+        vertex_to_neighbor_subdomain.clear();
+        vertex_to_neighbor_subdomain.resize(tria->n_vertices());
+        for (const auto &cell : tria->active_cell_iterators())
+          {
+            if (cell->is_ghost())
+              for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
+                vertex_to_neighbor_subdomain[cell->vertex_index(v)].insert(
+                  cell->subdomain_id());
+          }
+      }
+    return vertex_to_neighbor_subdomain;
+  }
+
 #include "grid_tools_cache.inst"
 
 } // namespace GridTools
index 83cfac4a13d997aa362828a3dd8e0f093467e489..3f3503aeadbac87f90e1fae0c5a8388de2064a98 100644 (file)
@@ -1186,16 +1186,8 @@ namespace Particles
         static_cast<typename std::vector<particle_iterator>::size_type>(
           particles.size() * 0.25));
 
-    std::vector<std::set<unsigned int>> vertex_to_neighbor_subdomain(
-      triangulation->n_vertices());
-
-    for (const auto &cell : triangulation->active_cell_iterators())
-      {
-        if (cell->is_ghost())
-          for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
-            vertex_to_neighbor_subdomain[cell->vertex_index(v)].insert(
-              cell->subdomain_id());
-      }
+    const std::vector<std::set<unsigned int>> vertex_to_neighbor_subdomain =
+      triangulation_cache->get_vertex_to_neighbor_subdomain();
 
     for (const auto &cell : triangulation->active_cell_iterators())
       {
@@ -1687,8 +1679,8 @@ namespace Particles
 
     switch (status)
       {
-        case parallel::distributed::Triangulation<dim, spacedim>::CELL_PERSIST:
-          {
+          case parallel::distributed::Triangulation<dim,
+                                                    spacedim>::CELL_PERSIST: {
             auto position_hint = particles.end();
             for (const auto &particle : loaded_particles_on_cell)
               {
@@ -1707,8 +1699,8 @@ namespace Particles
           }
           break;
 
-        case parallel::distributed::Triangulation<dim, spacedim>::CELL_COARSEN:
-          {
+          case parallel::distributed::Triangulation<dim,
+                                                    spacedim>::CELL_COARSEN: {
             typename std::multimap<internal::LevelInd,
                                    Particle<dim, spacedim>>::iterator
               position_hint = particles.end();
@@ -1733,8 +1725,8 @@ namespace Particles
           }
           break;
 
-        case parallel::distributed::Triangulation<dim, spacedim>::CELL_REFINE:
-          {
+          case parallel::distributed::Triangulation<dim,
+                                                    spacedim>::CELL_REFINE: {
             std::vector<
               typename std::multimap<internal::LevelInd,
                                      Particle<dim, spacedim>>::iterator>

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.