]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add another constructor for RemotePointEvaluation. 17126/head
authorDavid Wells <drwells@email.unc.edu>
Fri, 14 Jun 2024 20:01:54 +0000 (16:01 -0400)
committerDavid Wells <drwells@email.unc.edu>
Fri, 14 Jun 2024 20:01:54 +0000 (16:01 -0400)
Rather than setting up its own cache, we can just re-use the existing one. Also
re-use cached values.

include/deal.II/base/mpi_remote_point_evaluation.h
source/base/mpi_remote_point_evaluation.cc
tests/remote_point_evaluation/mapping_03.cc
tests/remote_point_evaluation/mapping_03.with_p4est=true.mpirun=2.output

index 9691d4349dc2f8d3fcddf2dfa22e921d2157b3f8..d7725d4dbc87059eeb450d41978093922395d799 100644 (file)
@@ -27,6 +27,9 @@ DEAL_II_NAMESPACE_OPEN
 // Forward declarations
 namespace GridTools
 {
+  template <int dim, int spacedim>
+  class Cache;
+
   namespace internal
   {
     template <int dim, int spacedim>
@@ -151,6 +154,20 @@ namespace Utilities
              const Triangulation<dim, spacedim> &tria,
              const Mapping<dim, spacedim>       &mapping);
 
+      /**
+       * Set up internal data structures and communication pattern based on
+       * an existing Cache and a list of points @p points.
+       *
+       * @warning This is a collective call that needs to be executed by all
+       *   processors in the communicator.
+       *
+       * @note If you want to be sure that all points have been found, call
+       *   all_points_found() after calling this function.
+       */
+      void
+      reinit(const GridTools::Cache<dim, spacedim> &cache,
+             const std::vector<Point<spacedim>>    &points);
+
       /**
        * Set up internal data structures and communication pattern based on
        * GridTools::internal::DistributedComputePointLocationsInternal.
index 76f6580d58a3c0896ca75dfbe1cefd035a506d4e..0789ae931449f69f074e0558b5804a262fb2e385 100644 (file)
@@ -88,32 +88,35 @@ namespace Utilities
       const Triangulation<dim, spacedim> &tria,
       const Mapping<dim, spacedim>       &mapping)
     {
+      const GridTools::Cache<dim, spacedim> cache(tria, mapping);
+
+      this->reinit(cache, points);
+    }
+
+
+
+    template <int dim, int spacedim>
+    void
+    RemotePointEvaluation<dim, spacedim>::reinit(
+      const GridTools::Cache<dim, spacedim> &cache,
+      const std::vector<Point<spacedim>>    &points)
+    {
 #ifndef DEAL_II_WITH_MPI
       Assert(false, ExcNeedsMPI());
+      (void)cache;
       (void)points;
-      (void)tria;
-      (void)mapping;
 #else
       if (tria_signal.connected())
         tria_signal.disconnect();
 
-      tria_signal =
-        tria.signals.any_change.connect([&]() { this->ready_flag = false; });
-
-      std::vector<BoundingBox<spacedim>> local_boxes;
-      for (const auto &cell :
-           tria.active_cell_iterators() | IteratorFilters::LocallyOwnedCell())
-        local_boxes.push_back(mapping.get_bounding_box(cell));
-
-      // create r-tree of bounding boxes
-      const auto local_tree = pack_rtree(local_boxes);
+      tria_signal = cache.get_triangulation().signals.any_change.connect(
+        [&]() { this->ready_flag = false; });
 
       // compress r-tree to a minimal set of bounding boxes
-      std::vector<std::vector<BoundingBox<spacedim>>> global_bboxes(1);
-      global_bboxes[0] =
-        extract_rtree_level(local_tree, additional_data.rtree_level);
-
-      const GridTools::Cache<dim, spacedim> cache(tria, mapping);
+      std::vector<std::vector<BoundingBox<spacedim>>> global_bboxes;
+      global_bboxes.emplace_back(
+        extract_rtree_level(cache.get_locally_owned_cell_bounding_boxes_rtree(),
+                            additional_data.rtree_level));
 
       const auto data =
         GridTools::internal::distributed_compute_point_locations(
@@ -126,7 +129,7 @@ namespace Utilities
           true,
           additional_data.enforce_unique_mapping);
 
-      this->reinit(data, tria, mapping);
+      this->reinit(data, cache.get_triangulation(), cache.get_mapping());
 #endif
     }
 
index 0b6364cdaf810ca63deedce7e8167948bd70a48c..86c1ec81c1c5786363689337f60511112fbc9f2d 100644 (file)
@@ -23,6 +23,7 @@
 #include <deal.II/fe/mapping_q1.h>
 
 #include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools_cache.h>
 
 #include <deal.II/numerics/vector_tools_evaluate.h>
 
@@ -32,7 +33,7 @@
 
 template <int dim>
 void
-test(const bool enforce_unique_map)
+test(const bool enforce_unique_map, const bool use_cache)
 {
   parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
   GridGenerator::subdivided_hyper_cube(tria, 2);
@@ -59,6 +60,11 @@ test(const bool enforce_unique_map)
   additional_data.tolerance              = 1e-6;
 
   Utilities::MPI::RemotePointEvaluation<dim> eval(additional_data);
+  if (use_cache)
+    {
+      GridTools::Cache<dim> cache(tria, mapping);
+      eval.reinit(cache, evaluation_points);
+    }
 
   const auto result_avg =
     VectorTools::point_values<1>(mapping,
@@ -89,6 +95,7 @@ main(int argc, char **argv)
   Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
   MPILogInitAll                    all;
 
-  test<2>(false);
-  test<2>(true);
+  test<2>(false, false);
+  test<2>(true, false);
+  test<2>(true, true);
 }
index 6fc910c6f3624d4ed8afe39476dbbfe43b80f4c6..1952162ffeb2772fe8de9398eaa09369cea0bc50 100644 (file)
@@ -27,6 +27,20 @@ DEAL:0::1e-11 1.00000 1.00000 1.00000 1.00000 1
 DEAL:0::1e-12 1.00000 1.00000 1.00000 1.00000 1
 DEAL:0::1e-13 1.00000 1.00000 1.00000 1.00000 1
 DEAL:0::1e-14 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::1e-1 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::1e-2 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::1e-3 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::1e-4 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::1e-5 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::1e-6 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::1e-7 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::1e-8 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::1e-9 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::1e-10 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::1e-11 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::1e-12 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::1e-13 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::1e-14 1.00000 1.00000 1.00000 1.00000 1
 
 DEAL:1::1e-1 2.00000 2.00000 2.00000 2.00000 2
 DEAL:1::1e-2 2.00000 2.00000 2.00000 2.00000 2
@@ -56,4 +70,18 @@ DEAL:1::1e-11 1.00000 1.00000 1.00000 1.00000 1
 DEAL:1::1e-12 1.00000 1.00000 1.00000 1.00000 1
 DEAL:1::1e-13 1.00000 1.00000 1.00000 1.00000 1
 DEAL:1::1e-14 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::1e-1 2.00000 2.00000 2.00000 2.00000 1
+DEAL:1::1e-2 2.00000 2.00000 2.00000 2.00000 1
+DEAL:1::1e-3 2.00000 2.00000 2.00000 2.00000 1
+DEAL:1::1e-4 2.00000 2.00000 2.00000 2.00000 1
+DEAL:1::1e-5 2.00000 2.00000 2.00000 2.00000 1
+DEAL:1::1e-6 2.00000 2.00000 2.00000 2.00000 1
+DEAL:1::1e-7 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::1e-8 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::1e-9 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::1e-10 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::1e-11 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::1e-12 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::1e-13 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::1e-14 1.00000 1.00000 1.00000 1.00000 1
 

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.