]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Adding global tree of bounding boxes 7001/head
authorGiovanni Alzetta <giovannialzetta@hotmail.it>
Mon, 30 Jul 2018 11:51:11 +0000 (13:51 +0200)
committerGiovanni Alzetta <giovannialzetta@hotmail.it>
Sat, 5 Jan 2019 09:39:14 +0000 (10:39 +0100)
doc/news/changes/minor/20181217GiovanniAlzetta [new file with mode: 0644]
include/deal.II/grid/grid_tools.h
source/grid/grid_tools.cc
source/grid/grid_tools.inst.in

diff --git a/doc/news/changes/minor/20181217GiovanniAlzetta b/doc/news/changes/minor/20181217GiovanniAlzetta
new file mode 100644 (file)
index 0000000..9e40685
--- /dev/null
@@ -0,0 +1,6 @@
+New: Function GridTools::build_global_description_tree which exchanges a given
+vector of bounding boxes on all processes and uses the result to build an Rtree
+with packing algorithm.
+<br>
+(Giovanni Alzetta, 2018/12/17)
+
index a9a1addc77490d0f996c96a99599392393a9299f..d6d87c9b759336896f45e606233bfba24ab22d38 100644 (file)
@@ -22,6 +22,8 @@
 #  include <deal.II/base/bounding_box.h>
 #  include <deal.II/base/geometry_info.h>
 
+#  include <deal.II/boost_adaptors/bounding_box.h>
+
 #  include <deal.II/dofs/dof_handler.h>
 
 #  include <deal.II/fe/mapping.h>
@@ -40,6 +42,7 @@
 
 #  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>
@@ -2690,6 +2693,36 @@ namespace GridTools
     const std::vector<BoundingBox<spacedim>> &local_bboxes,
     MPI_Comm                                  mpi_communicator);
 
+  /**
+   * In this collective operation each process provides a vector
+   * of bounding boxes and a communicator.
+   * All these vectors are gathered in each process
+   * and organized in a search tree which is then returned.
+   *
+   * The idea is that the vector of bounding boxes describes
+   * a relevant property which could be of use to other processes,
+   * e.g. for a distributed triangulation, the bounding
+   * boxes could describe the portion of the mesh which is
+   * locally owned by the current process.
+   *
+   * The search tree is an r-tree with packing algorithm,
+   * which is provided by boost library.
+   *
+   * In the returned tree, each node contains a pair of elements:
+   * the first being a bounding box,
+   * the second being the rank of the process whose local description
+   * contains the bounding box.
+   *
+   * Note: this function is a collective operation.
+   *
+   * @author Giovanni Alzetta, 2018.
+   */
+  template <int spacedim>
+  RTree<std::pair<BoundingBox<spacedim>, unsigned int>>
+  build_global_description_tree(
+    const std::vector<BoundingBox<spacedim>> &local_description,
+    MPI_Comm                                  mpi_communicator);
+
   /**
    * A structure that allows the transfer of cell data of type @p T from one processor
    * to another. It corresponds to a packed buffer that stores a vector of
index 14a856ec18d1025c051cee7c7ca10acb1cfebd49..2da178fa23efda55975c49680e902b86f495b7c0 100644 (file)
@@ -5081,6 +5081,65 @@ namespace GridTools
 #endif // DEAL_II_WITH_MPI
   }
 
+  template <int spacedim>
+  RTree<std::pair<BoundingBox<spacedim>, unsigned int>>
+  build_global_description_tree(
+    const std::vector<BoundingBox<spacedim>> &local_description,
+    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>>{};
+#else
+    // Exchanging local bounding boxes
+    std::vector<std::vector<BoundingBox<spacedim>>> global_bboxes;
+    global_bboxes =
+      Utilities::MPI::all_gather(mpi_communicator, local_description);
+
+    // Preparing to flatten the vector
+    unsigned int n_procs = Utilities::MPI::n_mpi_processes(mpi_communicator);
+    unsigned int tot_bboxes = 0;
+    // The i element of the following vector contains the index of the first
+    // local bounding box from the process of rank i
+    std::vector<unsigned int> bboxes_position(n_procs);
+
+    for (const auto &process_bboxes : global_bboxes)
+      tot_bboxes += process_bboxes.size();
+
+    // Now flattening the vector
+    std::vector<std::pair<BoundingBox<spacedim>, unsigned int>>
+      flat_global_bboxes;
+    flat_global_bboxes.reserve(tot_bboxes);
+    unsigned int process_index = 0;
+    for (const auto &process_bboxes : global_bboxes)
+      {
+        // Initializing a vector containing bounding boxes and rank of a process
+        std::vector<std::pair<BoundingBox<spacedim>, unsigned int>>
+          boxes_and_indices(process_bboxes.size());
+        // Adding to each box the rank of the process owning it
+        for (unsigned int i = 0; i < process_bboxes.size(); ++i)
+          boxes_and_indices[i] =
+            std::make_pair(process_bboxes[i], process_index);
+
+        flat_global_bboxes.insert(flat_global_bboxes.end(),
+                                  boxes_and_indices.begin(),
+                                  boxes_and_indices.end());
+        ++process_index;
+      }
+
+    // Building a tree out of the bounding boxes
+    // We avoid using the insert method so that boost uses the packing algorithm
+    RTree<std::pair<BoundingBox<spacedim>, unsigned int>> r_tree(
+      flat_global_bboxes.begin(), flat_global_bboxes.end());
+    return r_tree;
+#endif // DEAL_II_WITH_MPI
+  }
+
+
 
 } /* namespace GridTools */
 
index d8c1a41e5fcb4a5d2406e6083bc82def2d6f4ddd..c8b0ce9120e1095c4bcd418326caae361e2ab886 100644 (file)
@@ -149,6 +149,11 @@ for (deal_II_space_dimension : SPACE_DIMENSIONS)
     GridTools::guess_point_owner(
       const std::vector<std::vector<BoundingBox<deal_II_space_dimension>>> &,
       const std::vector<Point<deal_II_space_dimension>> &);
+
+    template RTree<
+      std::pair<BoundingBox<deal_II_space_dimension>, unsigned int>>
+    GridTools::build_global_description_tree(
+      const std::vector<BoundingBox<deal_II_space_dimension>> &, MPI_Comm);
   }
 
 
@@ -157,7 +162,6 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
 #if deal_II_dimension <= deal_II_space_dimension
     namespace GridTools
     \{
-
       template double
       diameter(
         const Triangulation<deal_II_dimension, deal_II_space_dimension> &);
@@ -331,9 +335,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
           &distorted_cells,
         Triangulation<deal_II_dimension, deal_II_space_dimension>
           &triangulation);
-
 #  endif
-
     \}
 #endif
   }

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.