]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added tests and support for non-mpi to GridTools::build_global_tree 7569/head
authorGiovanni Alzetta <giovannialzetta@hotmail.it>
Mon, 7 Jan 2019 11:09:28 +0000 (12:09 +0100)
committerGiovanni Alzetta <giovannialzetta@hotmail.it>
Wed, 9 Jan 2019 12:45:04 +0000 (13:45 +0100)
source/grid/grid_tools.cc
tests/grid/build_global_description_tree.cc [new file with mode: 0644]
tests/grid/build_global_description_tree.output [new file with mode: 0644]
tests/grid/build_global_description_tree.with_mpi=true.with_p4est=true.mpirun=2.output [new file with mode: 0644]

index 2da178fa23efda55975c49680e902b86f495b7c0..caa15ded21db24226353d16e62e0be617451f367 100644 (file)
@@ -5088,12 +5088,14 @@ namespace GridTools
     MPI_Comm                                  mpi_communicator)
   {
 #ifndef DEAL_II_WITH_MPI
-    (void)local_description;
     (void)mpi_communicator;
-    Assert(false,
-           ExcMessage(
-             "GridTools::build_global_description_tree() requires MPI."));
-    return RTree<std::pair<BoundingBox<spacedim>, unsigned int>>{};
+    // Building a tree with the only boxes available without MPI
+    std::vector<std::pair<BoundingBox<spacedim>, unsigned int>> boxes_index(
+      local_description.size());
+    // Adding to each box the rank of the process owning it
+    for (unsigned int i = 0; i < local_description.size(); ++i)
+      boxes_index[i] = std::make_pair(local_description[i], 0u);
+    return pack_rtree(boxes_index);
 #else
     // Exchanging local bounding boxes
     std::vector<std::vector<BoundingBox<spacedim>>> global_bboxes;
diff --git a/tests/grid/build_global_description_tree.cc b/tests/grid/build_global_description_tree.cc
new file mode 100644 (file)
index 0000000..3ece7f1
--- /dev/null
@@ -0,0 +1,124 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018-2019 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 GridTools::build_global_description_tree
+
+
+#include <deal.II/base/bounding_box.h>
+#include <deal.II/base/mpi.h>
+#include <deal.II/base/patterns.h>
+#include <deal.II/base/point.h>
+
+#include <deal.II/boost_adaptors/bounding_box.h>
+#include <deal.II/boost_adaptors/point.h>
+
+#include <deal.II/grid/grid_tools.h>
+
+#include <deal.II/numerics/rtree.h>
+
+#include <boost/archive/binary_iarchive.hpp>
+#include <boost/archive/binary_oarchive.hpp>
+#include <boost/geometry/index/detail/serialization.hpp>
+#include <boost/optional.hpp>
+#include <boost/serialization/array.hpp>
+#include <boost/serialization/vector.hpp>
+#include <boost/signals2.hpp>
+
+#include <algorithm>
+#include <utility>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+namespace bg  = boost::geometry;
+namespace bgi = boost::geometry::index;
+
+
+template <int dim, int spacedim>
+void
+test()
+{
+  unsigned int n_procs      = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+  unsigned int current_proc = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+  unsigned int next_p       = (current_proc + n_procs + 1) % n_procs;
+
+  // For this test we actually set up a simple system of bounding boxes
+  // The triangulation is really used only to get the right mpi_communicator
+  Point<spacedim> p1;
+  Point<spacedim> p2;
+  for (unsigned int d = 1; d < spacedim; ++d)
+    {
+      p1[d] = current_proc;
+      p2[d] = current_proc + 1;
+    }
+
+  BoundingBox<spacedim> box({p1, p2});
+
+  std::vector<BoundingBox<spacedim>> local_description(1, box);
+
+  auto built_tree =
+    GridTools::build_global_description_tree(local_description, MPI_COMM_WORLD);
+
+  Point<spacedim> my_point;
+  Point<spacedim> point_inside_d_1;
+  Point<spacedim> outside_point;
+
+  for (unsigned int d = 1; d < spacedim; ++d)
+    {
+      my_point[d]         = current_proc + 0.1 + d / 4.2;
+      point_inside_d_1[d] = next_p + 0.25 + d / 5.1;
+      outside_point[d]    = -current_proc - 2.6;
+    }
+
+  std::vector<std::pair<BoundingBox<spacedim>, unsigned int>> test_results;
+  built_tree.query(bgi::intersects(outside_point),
+                   std::back_inserter(test_results));
+  if (test_results.size() != 0)
+    deallog
+      << "Point found inside a bounding box! It should be outside all of them!"
+      << std::endl;
+
+  built_tree.query(bgi::intersects(my_point), std::back_inserter(test_results));
+  if (std::get<1>(test_results[0]) != current_proc)
+    deallog << "Error: Point found inside wrong process: "
+            << std::get<1>(test_results[0]) << std::endl;
+
+  built_tree.query(bgi::intersects(point_inside_d_1),
+                   std::back_inserter(test_results));
+  if (std::get<1>(test_results[1]) != next_p)
+    deallog << "Error: Point found inside wrong process: "
+            << std::get<1>(test_results[1]) << std::endl;
+}
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+#ifdef DEAL_II_WITH_MPI
+  MPILogInitAll log;
+#else
+  initlog();
+  deallog.push("0");
+#endif
+
+  deallog << "Deal.II build_global_description_tree:" << std::endl;
+
+  test<2, 2>();
+  test<3, 3>();
+
+  deallog << "Test finished" << std::endl;
+  return 0;
+}
diff --git a/tests/grid/build_global_description_tree.output b/tests/grid/build_global_description_tree.output
new file mode 100644 (file)
index 0000000..7958569
--- /dev/null
@@ -0,0 +1,3 @@
+
+DEAL:0::Deal.II build_global_description_tree:
+DEAL:0::Test finished
diff --git a/tests/grid/build_global_description_tree.with_mpi=true.with_p4est=true.mpirun=2.output b/tests/grid/build_global_description_tree.with_mpi=true.with_p4est=true.mpirun=2.output
new file mode 100644 (file)
index 0000000..b16cf2d
--- /dev/null
@@ -0,0 +1,7 @@
+
+DEAL:0::Deal.II build_global_description_tree:
+DEAL:0::Test finished
+
+DEAL:1::Deal.II build_global_description_tree:
+DEAL:1::Test finished
+

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.