]> https://gitweb.dealii.org/ - dealii.git/commitdiff
RPE: refactor determination of unique mapping and all points found 13604/head
authorPeter Munch <peterrmuench@gmail.com>
Sat, 9 Apr 2022 09:44:03 +0000 (11:44 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Wed, 13 Apr 2022 07:29:21 +0000 (09:29 +0200)
include/deal.II/base/mpi_remote_point_evaluation.h
source/base/mpi_remote_point_evaluation.cc
tests/remote_point_evaluation/mapping_04.cc [new file with mode: 0644]
tests/remote_point_evaluation/mapping_04.with_p4est=true.mpirun=2.output [new file with mode: 0644]

index cf1429a536268db6baba284ba9dc821ded312ac4..51a0a757782dc0352abeacb44375c2d9697e8efb 100644 (file)
@@ -242,6 +242,11 @@ namespace Utilities
        */
       bool unique_mapping;
 
+      /**
+       * Cache if all points passed in during reinit() have been found.
+       */
+      bool all_points_found_flag;
+
       /**
        * Since for each point multiple or no results can be available, the
        * pointers in this vector indicate the first and last entry associated
index dd1907609a9493240168c099f9d82466c64ec4a8..671cd9e4671a5aae8ad6e737f3d6480419b6c539 100644 (file)
@@ -130,14 +130,52 @@ namespace Utilities
           this->point_ptrs[std::get<1>(data.recv_components[i]) + 1]++;
         }
 
-      unique_mapping = true;
+      std::tuple<unsigned int, unsigned int> n_owning_processes_default{
+        numbers::invalid_unsigned_int, 0};
+      std::tuple<unsigned int, unsigned int> n_owning_processes_local =
+        n_owning_processes_default;
+
       for (unsigned int i = 0; i < points.size(); ++i)
         {
-          if (unique_mapping && this->point_ptrs[i + 1] != 1)
-            unique_mapping = false;
+          std::get<0>(n_owning_processes_local) =
+            std::min(std::get<0>(n_owning_processes_local),
+                     this->point_ptrs[i + 1]);
+          std::get<1>(n_owning_processes_local) =
+            std::max(std::get<1>(n_owning_processes_local),
+                     this->point_ptrs[i + 1]);
+
           this->point_ptrs[i + 1] += this->point_ptrs[i];
         }
 
+      const auto n_owning_processes_global =
+        Utilities::MPI::all_reduce<std::tuple<unsigned int, unsigned int>>(
+          n_owning_processes_local,
+          tria.get_communicator(),
+          [&](const auto &a,
+              const auto &b) -> std::tuple<unsigned int, unsigned int> {
+            if (a == n_owning_processes_default)
+              return b;
+
+            if (b == n_owning_processes_default)
+              return a;
+
+            return std::tuple<unsigned int, unsigned int>{
+              std::min(std::get<0>(a), std::get<0>(b)),
+              std::max(std::get<1>(a), std::get<1>(b))};
+          });
+
+      if (n_owning_processes_global == n_owning_processes_default)
+        {
+          unique_mapping        = true;
+          all_points_found_flag = true;
+        }
+      else
+        {
+          unique_mapping = (std::get<0>(n_owning_processes_global) == 1) &&
+                           (std::get<1>(n_owning_processes_global) == 1);
+          all_points_found_flag = std::get<0>(n_owning_processes_global) > 0;
+        }
+
       Assert(enforce_unique_mapping == false || unique_mapping,
              ExcInternalError());
 
@@ -189,15 +227,7 @@ namespace Utilities
     bool
     RemotePointEvaluation<dim, spacedim>::all_points_found() const
     {
-      if (is_map_unique())
-        return true;
-
-      if (point_ptrs.size() > 0)
-        for (unsigned int i = 0; i < point_ptrs.size() - 1; ++i)
-          if (point_found(i) == false)
-            return false;
-
-      return true;
+      return all_points_found_flag;
     }
 
 
@@ -208,7 +238,11 @@ namespace Utilities
       const unsigned int i) const
     {
       AssertIndexRange(i, point_ptrs.size() - 1);
-      return (point_ptrs[i + 1] - point_ptrs[i]) > 0;
+
+      if (all_points_found_flag)
+        return true;
+      else
+        return (point_ptrs[i + 1] - point_ptrs[i]) > 0;
     }
 
 
diff --git a/tests/remote_point_evaluation/mapping_04.cc b/tests/remote_point_evaluation/mapping_04.cc
new file mode 100644 (file)
index 0000000..99c9d26
--- /dev/null
@@ -0,0 +1,197 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// Test Utilities::MPI::RemotePointEvaluation if points have been found
+// and if the mapping is unique.
+
+#include <deal.II/base/mpi.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/mapping_q1.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/lac/la_vector.h>
+
+#include <deal.II/numerics/vector_tools_evaluate.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+do_test(const bool                     enforce_unique_map,
+        const std::vector<Point<dim>> &evaluation_points,
+        const std::pair<bool, bool> &  expected_result)
+{
+  parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+  GridGenerator::subdivided_hyper_rectangle(tria, {2, 1}, {0, 0}, {2, 1});
+
+  MappingQ1<dim> mapping;
+
+  Utilities::MPI::RemotePointEvaluation<dim> eval(1.e-6, enforce_unique_map);
+
+  eval.reinit(evaluation_points, tria, mapping);
+
+  Assert(eval.is_map_unique() == expected_result.first, ExcInternalError());
+  Assert(eval.all_points_found() == expected_result.second, ExcInternalError());
+}
+
+
+template <int dim>
+void
+test()
+{
+  const auto my_rank = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+
+  {
+    // (-, -)
+    std::vector<Point<dim>> evaluation_points;
+
+    do_test(false, evaluation_points, {true, true});
+  }
+
+  {
+    // (not found, -)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 0)
+      evaluation_points.emplace_back(2.5, 0.5);
+
+    do_test(false, evaluation_points, {false, false});
+  }
+
+  {
+    // (unique, -)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 0)
+      evaluation_points.emplace_back(1.5, 0.5);
+
+    do_test(false, evaluation_points, {true, true});
+  }
+
+  {
+    // (non-unique, -)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 0)
+      evaluation_points.emplace_back(1.0, 0.5);
+
+    do_test(false, evaluation_points, {false, true});
+  }
+
+  {
+    // (enforced-unique, -)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 0)
+      evaluation_points.emplace_back(1.0, 0.5);
+
+    do_test(true, evaluation_points, {true, true});
+  }
+
+  {
+    // (-, not found)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 1)
+      evaluation_points.emplace_back(2.5, 0.5);
+
+    do_test(false, evaluation_points, {false, false});
+  }
+
+  {
+    // (-, unique)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 1)
+      evaluation_points.emplace_back(0.5, 0.5);
+
+    do_test(false, evaluation_points, {true, true});
+  }
+
+  {
+    // (-, non-unique)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 1)
+      evaluation_points.emplace_back(1.0, 0.5);
+
+    do_test(false, evaluation_points, {false, true});
+  }
+
+  {
+    // (-, enforced-unique)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 1)
+      evaluation_points.emplace_back(1.0, 0.5);
+
+    do_test(true, evaluation_points, {true, true});
+  }
+
+  {
+    // (unique, not found)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 0)
+      evaluation_points.emplace_back(1.5, 0.5);
+    if (my_rank == 1)
+      evaluation_points.emplace_back(2.5, 0.5);
+
+    do_test(false, evaluation_points, {false, false});
+  }
+
+  {
+    // (unique, not unique)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 0)
+      evaluation_points.emplace_back(1.5, 0.5);
+    if (my_rank == 1)
+      evaluation_points.emplace_back(1.0, 0.5);
+
+    do_test(false, evaluation_points, {false, true});
+  }
+
+  {
+    // (unique, enforced unique)
+    std::vector<Point<dim>> evaluation_points;
+
+    if (my_rank == 0)
+      evaluation_points.emplace_back(1.5, 0.5);
+    if (my_rank == 1)
+      evaluation_points.emplace_back(1.0, 0.5);
+
+    do_test(true, evaluation_points, {true, true});
+  }
+}
+
+
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  AssertDimension(Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD), 2);
+
+  test<2>();
+}
diff --git a/tests/remote_point_evaluation/mapping_04.with_p4est=true.mpirun=2.output b/tests/remote_point_evaluation/mapping_04.with_p4est=true.mpirun=2.output
new file mode 100644 (file)
index 0000000..b28b04f
--- /dev/null
@@ -0,0 +1,3 @@
+
+
+

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.