]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Completely reworked the pull request.
authorBruno Blais <blais.bruno@gmail.com>
Thu, 12 Mar 2020 13:03:18 +0000 (09:03 -0400)
committerBruno Blais <blais.bruno@gmail.com>
Thu, 12 Mar 2020 13:03:18 +0000 (09:03 -0400)
The particle handler now has a GridTools::Cache that is used to obtain the vertex_to_cell_centers and vertex_to_cells
This enables the same feature as before but is a lot cleaner and flexible.
However, the cache has to be stored in a unique pointer because the particlehandler has a constructor that enables it to be constructed without a triangulation.
This is not the case for the GridTools::Cache. Consequently, I thought that unique an std::unique was a good compromise here.

As always, looking forward to your comments :)!

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

index 5fec6138200054830d6c37ba02a8c020bc68c0a3..78f9d22f504581a079dfec5577dc5166b71859fd 100644 (file)
@@ -29,6 +29,8 @@
 
 #include <deal.II/fe/mapping.h>
 
+#include <deal.II/grid/grid_tools_cache.h>
+
 #include <deal.II/particles/particle.h>
 #include <deal.II/particles/particle_iterator.h>
 #include <deal.II/particles/property_pool.h>
@@ -86,15 +88,12 @@ namespace Particles
     /**
      * Destructor.
      */
-    virtual ~ParticleHandler() override
-    {
-      connection.disconnect();
-    }
+    virtual ~ParticleHandler() override = default;
 
     /**
      * Initialize the particle handler. This function does not clear the
-     * internal data structures, it just sets the triangulation and the mapping
-     * to be used.
+     * internal data structures, it just sets the triangulation and the
+     * mapping to be used.
      */
     void
     initialize(const Triangulation<dim, spacedim> &tria,
@@ -598,16 +597,6 @@ namespace Particles
     void
     serialize(Archive &ar, const unsigned int version);
 
-    /**
-     * Updates the structures that were generated to search in which cells the
-     * particles belong to. The call to this function updates the
-     * vertex_to_cells and the vertex_to_cell_centers data structure. This
-     * function is called every time the triangulation is altered.
-     *
-     */
-    void
-    update_triangulation_cache();
-
   private:
     /**
      * Address of the triangulation to work on.
@@ -709,30 +698,15 @@ namespace Particles
      */
     unsigned int handle;
 
-    /**
-     * A connection to the corresponding any_changes signal of the Triangulation
-     * which is attached to the particle_handler
-     */
-    boost::signals2::connection connection;
-
-    /** This variable stores a set from vertices to adjacent cells. It is used
-     * to locate the cells in which the particle reside. This structured is kept
-     * in memory to prevent its recalculation every time we
-     * sort_into_subdomain_and_cells(). This structure is updated everytime
-     * we call update_triangulation_cache.
-     */
-    std::vector<
-      std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
-      vertex_to_cells;
-
-
-    /** This variable stores a vector of vectors to cell centers. It is used
-     * to locate the cells in which the particle reside. This structured is kept
-     * in memory to prevent its recalculation every time we
-     * sort_into_subdomain_and_cells(). This structure is updated everytime
-     * we call update_triangulation_cache.
+    /** The GridTools::Cache is used to store the information about the
+     * vertex_to_cells set and the vertex_to_cell_centers vectors to prevent
+     * recomputing them every time we sort_into_subdomain_and_cells()
+     * This cache is automatically updated when the triangulation has
+     * changed. This cache is stored within a unique pointer because the
+     * particle handler has a constructor that enables it to be constructed
+     * without a triangulation. The cache does not have such a constructor.
      */
-    std::vector<std::vector<Tensor<1, spacedim>>> vertex_to_cell_centers;
+    std::unique_ptr<GridTools::Cache<dim, spacedim>> triangulation_cache;
 
 #ifdef DEAL_II_WITH_MPI
     /**
index 90b8729b504c01754789ec21297c6ebfbae82a27..880b85b185e14332eb14a1df09a9ac4e5fbc276e 100644 (file)
@@ -120,13 +120,9 @@ namespace Particles
     , load_callback()
     , handle(numbers::invalid_unsigned_int)
   {
-    // Update the triangulation cache to ensure that you can sort the particles
-    // correctly
-    update_triangulation_cache();
-
-    connection = triangulation.signals.any_change.connect(
-      std::bind(&ParticleHandler<dim, spacedim>::update_triangulation_cache,
-                std::ref(*this)));
+    triangulation_cache =
+      std_cxx14::make_unique<GridTools::Cache<dim, spacedim>>(triangulation,
+                                                              mapping);
   }
 
 
@@ -144,13 +140,11 @@ namespace Particles
     // Create the memory pool that will store all particle properties
     property_pool = std_cxx14::make_unique<PropertyPool>(n_properties);
 
-    // Update the triangulation cache to ensure that you can sort the particles
-    // correctly
-    update_triangulation_cache();
-
-    connection = triangulation->signals.any_change.connect(
-      std::bind(&ParticleHandler<dim, spacedim>::update_triangulation_cache,
-                std::ref(*this)));
+    // Create the grid cache to cache the informations about the triangulation
+    // that is used to locate the particles into subdomains and cells
+    triangulation_cache =
+      std_cxx14::make_unique<GridTools::Cache<dim, spacedim>>(new_triangulation,
+                                                              new_mapping);
   }
 
 
@@ -926,7 +920,16 @@ namespace Particles
         static_cast<vector_size>(particles_out_of_cell.size() * 0.25));
 
     {
-      // update_triangulation_cache();
+      // Create a map from vertices to adjacent cells using grid cache
+      std::vector<
+        std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
+        vertex_to_cells = triangulation_cache->get_vertex_to_cell_map();
+
+      // Create a corresponding map of vectors from vertex to cell center using
+      // grid cache
+      std::vector<std::vector<Tensor<1, spacedim>>> vertex_to_cell_centers =
+        triangulation_cache->get_vertex_to_cell_centers_directions();
+
       std::vector<unsigned int> neighbor_permutation;
 
       // Find the cells that the particles moved to.
@@ -1751,21 +1754,6 @@ namespace Particles
           break;
       }
   }
-
-  template <int dim, int spacedim>
-  void
-  ParticleHandler<dim, spacedim>::update_triangulation_cache()
-  {
-    // Create a map from vertices to adjacent cells
-    vertex_to_cells = GridTools::vertex_to_cell_map(*triangulation);
-
-    // Create a corresponding map of vectors from vertex to cell center
-    vertex_to_cell_centers =
-      GridTools::vertex_to_cell_centers_directions(*triangulation,
-                                                   vertex_to_cells);
-  }
-
-
 } // namespace Particles
 
 #include "particle_handler.inst"

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.