]> https://gitweb.dealii.org/ - dealii.git/commitdiff
First version of an update_ghost_particles
authorBruno <blais.bruno@gmail.com>
Sat, 31 Oct 2020 16:25:07 +0000 (12:25 -0400)
committerBruno <blais.bruno@gmail.com>
Sat, 31 Oct 2020 16:25:07 +0000 (12:25 -0400)
Does not contain any modification to send_recv

include/deal.II/particles/particle_handler.h
source/particles/particle_handler.cc

index 5752de864e5c6d53b838d4cdee1e7883370161ff..897adad9970f04e7ed7e62e9d88413eff1fefcbd 100644 (file)
@@ -707,6 +707,13 @@ namespace Particles
     void
     exchange_ghost_particles();
 
+    /**
+     * Update all particles that live in cells that are ghost cells to
+     * other processes.
+     */
+    void
+    update_ghost_particles();
+
     /**
      * Callback function that should be called before every refinement
      * and when writing checkpoints. This function is used to
@@ -809,6 +816,14 @@ namespace Particles
      */
     std::multimap<internal::LevelInd, Particle<dim, spacedim>> ghost_particles;
 
+    /**
+     * Set of particles that currently live in the ghost cells of the local
+     * domain, organized by the subdomain_id. These
+     * particles are equivalent to the ghost entries in distributed vectors.
+     */
+    std::map<types::subdomain_id, std::vector<particle_iterator>>
+      ghost_particles_by_domain;
+
     /**
      * This variable stores how many particles are stored globally. It is
      * calculated by update_cached_numbers().
index a41ad79c2afbbd0a63cb7a9589ac50a78b8a55aa..f0e665db428185a84c6d889714508f57ef8c1ddf 100644 (file)
@@ -90,6 +90,7 @@ namespace Particles
     , property_pool(std::make_unique<PropertyPool>(0))
     , particles()
     , ghost_particles()
+    , ghost_particles_by_domain()
     , global_number_of_particles(0)
     , global_max_particles_per_cell(0)
     , next_free_particle_index(0)
@@ -111,6 +112,7 @@ namespace Particles
     , property_pool(std::make_unique<PropertyPool>(n_properties))
     , particles()
     , ghost_particles()
+    , ghost_particles_by_domain()
     , global_number_of_particles(0)
     , global_max_particles_per_cell(0)
     , next_free_particle_index(0)
@@ -164,10 +166,11 @@ namespace Particles
     global_number_of_particles = particle_handler.global_number_of_particles;
     global_max_particles_per_cell =
       particle_handler.global_max_particles_per_cell;
-    next_free_particle_index = particle_handler.next_free_particle_index;
-    particles                = particle_handler.particles;
-    ghost_particles          = particle_handler.ghost_particles;
-    handle                   = particle_handler.handle;
+    next_free_particle_index  = particle_handler.next_free_particle_index;
+    particles                 = particle_handler.particles;
+    ghost_particles           = particle_handler.ghost_particles;
+    ghost_particles_by_domain = particle_handler.ghost_particles_by_domain;
+    handle                    = particle_handler.handle;
 
     // copy dynamic properties
     auto from_particle = particle_handler.begin();
@@ -1176,8 +1179,8 @@ namespace Particles
     // First clear the current ghost_particle information
     ghost_particles.clear();
 
-    std::map<types::subdomain_id, std::vector<particle_iterator>>
-      ghost_particles_by_domain;
+    ghost_particles_by_domain.clear();
+
 
     const std::set<types::subdomain_id> ghost_owners =
       parallel_triangulation->ghost_owners();
@@ -1222,6 +1225,31 @@ namespace Particles
 #endif
   }
 
+  template <int dim, int spacedim>
+  void
+  ParticleHandler<dim, spacedim>::update_ghost_particles()
+  {
+    // Nothing to do in serial computations
+    const auto parallel_triangulation =
+      dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
+        &*triangulation);
+    if (parallel_triangulation != nullptr)
+      {
+        if (dealii::Utilities::MPI::n_mpi_processes(
+              parallel_triangulation->get_communicator()) == 1)
+          return;
+      }
+    else
+      return;
+
+#ifdef DEAL_II_WITH_MPI
+    // First clear the current ghost_particle information
+    ghost_particles.clear();
+
+    send_recv_particles(ghost_particles_by_domain, ghost_particles);
+#endif
+  }
+
 
 
 #ifdef DEAL_II_WITH_MPI
@@ -1679,8 +1707,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)
               {
@@ -1699,8 +1727,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();
@@ -1725,8 +1753,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.