]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add AdditionalData struct to RemotePointEvaluation 16323/head
authorMagdalena Schreter <magdalena.schreter@tum.de>
Wed, 6 Dec 2023 15:20:19 +0000 (16:20 +0100)
committerMagdalena Schreter <magdalena.schreter@tum.de>
Sun, 31 Dec 2023 13:13:36 +0000 (14:13 +0100)
doc/news/changes/minor/20231206Schreter [new file with mode: 0644]
include/deal.II/base/mpi_remote_point_evaluation.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
source/base/mpi_remote_point_evaluation.cc
tests/remote_point_evaluation/mapping_01.cc
tests/remote_point_evaluation/mapping_03.cc
tests/remote_point_evaluation/mapping_04.cc
tests/remote_point_evaluation/search_adjacent_cells.cc

diff --git a/doc/news/changes/minor/20231206Schreter b/doc/news/changes/minor/20231206Schreter
new file mode 100644 (file)
index 0000000..f3c7df1
--- /dev/null
@@ -0,0 +1,6 @@
+New: Add AdditionalData struct for configuring the RemotePointEvaluation class. 
+The latter is passed into a new constructor.
+Deprecated: The old constructor of RemotePointEvaluation class taking the 
+parameters value-by-value is marked as deprecated.
+<br>
+(Magdalena Schreter-Fleischhacker, Peter Munch, 2023/12/06)
index 0c441c2854a2b5b92b210c07ac7f127a9ed4bbde..8f9934e2d3b0556f7b4ef7bb0a7ffca18d0298e9 100644 (file)
@@ -53,9 +53,62 @@ namespace Utilities
     class RemotePointEvaluation
     {
     public:
+      /**
+       * AdditionalData structure that can be used to tweak parameters
+       * of RemotePointEvaluation.
+       */
+      struct AdditionalData
+      {
+      public:
+        /**
+         *  Constructor.
+         */
+        AdditionalData(
+          const double       tolerance                              = 1e-6,
+          const bool         enforce_unique_mapping                 = false,
+          const unsigned int rtree_level                            = 0,
+          const std::function<std::vector<bool>()> &marked_vertices = {});
+
+        /**
+         * Tolerance in terms of unit cell coordinates for determining all cells
+         * around a point passed to RemotePointEvaluation during reinit().
+         * Depending on the problem, it might be necessary to adjust the
+         * tolerance in order to be able to identify a cell. Floating point
+         * arithmetic implies that a point will, in general, not lie exactly on
+         * a vertex, edge, or face.
+         */
+        double tolerance;
+
+        /**
+         * Enforce unique mapping, i.e., (one-to-one) relation of points and
+         * cells.
+         */
+        bool enforce_unique_mapping;
+
+        /**
+         * RTree level to be used during the construction of the bounding boxes.
+         */
+        unsigned int rtree_level;
+
+        /**
+         * Function that marks relevant vertices to make search of active cells
+         * around point more efficient.
+         */
+        std::function<std::vector<bool>()> marked_vertices;
+      };
+
       /**
        * Constructor.
        *
+       * @param additional_data Configure options for RemotePointEvaluation.
+       */
+      RemotePointEvaluation(
+        const AdditionalData &additional_data = AdditionalData());
+
+      /**
+       * Constructor. This constructor is deprecated. Use the other constructor
+       * taking AdditionalData instead.
+       *
        * @param tolerance Tolerance in terms of unit cell coordinates for
        *   determining all cells around a point passed to the class during
        *   reinit(). Depending on the problem, it might be necessary to adjust
@@ -67,9 +120,13 @@ namespace Utilities
        * @param rtree_level RTree level to be used during the construction of the bounding boxes.
        * @param marked_vertices Function that marks relevant vertices to make search
        *   of active cells around point more efficient.
+       *
+       * @deprecated
        */
+      DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
+        "Use the constructor with AdditionalData struct.")
       RemotePointEvaluation(
-        const double       tolerance                              = 1e-6,
+        const double       tolerance,
         const bool         enforce_unique_mapping                 = false,
         const unsigned int rtree_level                            = 0,
         const std::function<std::vector<bool>()> &marked_vertices = {});
@@ -312,27 +369,9 @@ namespace Utilities
 
     private:
       /**
-       * Tolerance to be used while determining the surrounding cells of a
-       * point.
-       */
-      const double tolerance;
-
-      /**
-       * Enforce unique mapping, i.e., (one-to-one) relation of points and
-       * cells.
-       */
-      const bool enforce_unique_mapping;
-
-      /**
-       * RTree level to be used during the construction of the bounding boxes.
-       */
-      const unsigned int rtree_level;
-
-      /**
-       * Function that marks relevant vertices to make search of active cells
-       * around point more efficient.
+       * Additional data with basic settings.
        */
-      const std::function<std::vector<bool>()> marked_vertices;
+      const AdditionalData additional_data;
 
       /**
        * Storage for the status of the triangulation signal.
index 65680d32e887c0d62cd10a0d1ded958bcb6af1aa..323e9c27cb8ee1b86a5ad242012b32dc389d0413 100644 (file)
@@ -4991,7 +4991,11 @@ template <int dim, typename Number>
 MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
   MGTwoLevelTransferNonNested(const AdditionalData &data)
   : additional_data(data)
-  , rpe(data.tolerance, false, data.rtree_level, {})
+  , rpe(typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData(
+      data.tolerance,
+      false,
+      data.rtree_level,
+      {}))
 {}
 
 template <int dim, typename Number>
index e0c52cceffd30f3caf72fd83298bba4e6fab4d4d..27232e1bdfd55c330a3c2cc4eee7fdf3c151d725 100644 (file)
@@ -36,7 +36,7 @@ namespace Utilities
   namespace MPI
   {
     template <int dim, int spacedim>
-    RemotePointEvaluation<dim, spacedim>::RemotePointEvaluation(
+    RemotePointEvaluation<dim, spacedim>::AdditionalData::AdditionalData(
       const double                              tolerance,
       const bool                                enforce_unique_mapping,
       const unsigned int                        rtree_level,
@@ -45,6 +45,29 @@ namespace Utilities
       , enforce_unique_mapping(enforce_unique_mapping)
       , rtree_level(rtree_level)
       , marked_vertices(marked_vertices)
+    {}
+
+
+
+    template <int dim, int spacedim>
+    RemotePointEvaluation<dim, spacedim>::RemotePointEvaluation(
+      const AdditionalData &additional_data)
+      : additional_data(additional_data)
+      , ready_flag(false)
+    {}
+
+
+
+    template <int dim, int spacedim>
+    RemotePointEvaluation<dim, spacedim>::RemotePointEvaluation(
+      const double                              tolerance,
+      const bool                                enforce_unique_mapping,
+      const unsigned int                        rtree_level,
+      const std::function<std::vector<bool>()> &marked_vertices)
+      : additional_data(tolerance,
+                        enforce_unique_mapping,
+                        rtree_level,
+                        marked_vertices)
       , ready_flag(false)
     {}
 
@@ -88,7 +111,8 @@ namespace Utilities
 
       // 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, rtree_level);
+      global_bboxes[0] =
+        extract_rtree_level(local_tree, additional_data.rtree_level);
 
       const GridTools::Cache<dim, spacedim> cache(tria, mapping);
 
@@ -97,10 +121,11 @@ namespace Utilities
           cache,
           points,
           global_bboxes,
-          marked_vertices ? marked_vertices() : std::vector<bool>(),
-          tolerance,
+          additional_data.marked_vertices ? additional_data.marked_vertices() :
+                                            std::vector<bool>(),
+          additional_data.tolerance,
           true,
-          enforce_unique_mapping);
+          additional_data.enforce_unique_mapping);
 
       this->reinit(data, tria, mapping);
 #endif
@@ -185,7 +210,7 @@ namespace Utilities
           all_points_found_flag = std::get<0>(n_owning_processes_global) > 0;
         }
 
-      Assert(enforce_unique_mapping == false || unique_mapping,
+      Assert(additional_data.enforce_unique_mapping == false || unique_mapping,
              ExcInternalError());
 
       cell_data        = std::make_unique<CellData>(tria);
index 9a650556798fab8c74b8423f770374270903e813..454c7fa4ff997072adaaa2b79b3675afd263d32c 100644 (file)
@@ -58,7 +58,11 @@ test(const bool enforce_unique_map)
     for (unsigned int i = 0; i <= 4; ++i)
       evaluation_points.emplace_back(i * 0.25, j * 0.25);
 
-  Utilities::MPI::RemotePointEvaluation<dim> eval(1e-6, enforce_unique_map);
+  typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+    additional_data;
+  additional_data.enforce_unique_mapping = enforce_unique_map;
+  additional_data.tolerance              = 1e-6;
+  Utilities::MPI::RemotePointEvaluation<dim> eval(additional_data);
 
   const auto result_avg =
     VectorTools::point_values<1>(mapping,
index 589530a332db0dbd160977aa2284df2206b5fe86..50a553d614a34aa221029c16d61567476bb5ea9e 100644 (file)
@@ -55,7 +55,12 @@ test(const bool enforce_unique_map)
     // should be assigned to rank 0 for unique mapping
     evaluation_points.emplace_back(0.5, 0.5 + std::pow<double>(10.0, -i));
 
-  Utilities::MPI::RemotePointEvaluation<dim> eval(1e-6, enforce_unique_map);
+  typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+    additional_data;
+  additional_data.enforce_unique_mapping = enforce_unique_map;
+  additional_data.tolerance              = 1e-6;
+
+  Utilities::MPI::RemotePointEvaluation<dim> eval(additional_data);
 
   const auto result_avg =
     VectorTools::point_values<1>(mapping,
index c791aa78dcb1febc59d76fe79d5f23edee8f559b..df321285a41a50ade571a579333bed1374003515 100644 (file)
@@ -42,7 +42,11 @@ do_test(const bool                     enforce_unique_map,
 
   MappingQ1<dim> mapping;
 
-  Utilities::MPI::RemotePointEvaluation<dim> eval(1.e-6, enforce_unique_map);
+  typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+    additional_data;
+  additional_data.enforce_unique_mapping = enforce_unique_map;
+
+  Utilities::MPI::RemotePointEvaluation<dim> eval(additional_data);
 
   eval.reinit(evaluation_points, tria, mapping);
 
index 8e95d38684a31184e7c16098fc75809b205c3338..458fe1b7baca9cc608974b9ad5481bd31558c333 100644 (file)
@@ -158,7 +158,10 @@ test(const unsigned int mapping_degree,
     }
 
   // initialize RPE without any marked points
-  dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe(tolerance);
+  typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+    additional_data;
+  additional_data.tolerance = tolerance;
+  dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe(additional_data);
   rpe.reinit(points, tria, mapping);
 
   unsigned int                    n_points_not_found_rpe = 0;
@@ -182,8 +185,13 @@ test(const unsigned int mapping_degree,
 
   // initialize RPE with all points marked
   std::vector<bool> marked_vertices(tria.n_vertices(), true);
-  dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe2(
-    tolerance, false, 0, [marked_vertices]() { return marked_vertices; });
+
+  typename Utilities::MPI::RemotePointEvaluation<dim>::AdditionalData
+    additional_data2(tolerance, false, 0, [marked_vertices]() {
+      return marked_vertices;
+    });
+
+  dealii::Utilities::MPI::RemotePointEvaluation<dim> rpe2(additional_data2);
 
   rpe2.reinit(points, tria, mapping);
 

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.