]> https://gitweb.dealii.org/ - dealii.git/commitdiff
New extract_rtree_level function, and RTree visitor.
authorLuca Heltai <luca.heltai@sissa.it>
Mon, 20 Jan 2020 14:26:42 +0000 (15:26 +0100)
committerLuca Heltai <luca.heltai@sissa.it>
Sat, 9 May 2020 22:49:30 +0000 (00:49 +0200)
doc/news/changes/minor/20200503LucaHeltai [new file with mode: 0644]
include/deal.II/base/patterns.h
include/deal.II/numerics/rtree.h
tests/boost/extract_reprsentative_set_01.cc [new file with mode: 0644]
tests/boost/extract_reprsentative_set_01.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20200503LucaHeltai b/doc/news/changes/minor/20200503LucaHeltai
new file mode 100644 (file)
index 0000000..a0379e6
--- /dev/null
@@ -0,0 +1,4 @@
+New: ExtractLevelVisitor and extract_rtree_level() allow one to return the vector of BoundingBox objects 
+that make up a specific level of a boost::rtree object.
+<br>
+(Luca Heltai, 2020/05/03)
index e37d5bb8fa4845f70c0c89b441d481f1bc788793..9541bbbf79aaceb510a1cc9e100dd2c253353ee4 100644 (file)
@@ -2231,15 +2231,10 @@ namespace Patterns
       {
         static_assert(internal::RankInfo<T>::map_rank > 0,
                       "Cannot use this class for non Map-compatible types.");
-        return std_cxx14::make_unique<Patterns::Map>(
+        return std_cxx14::make_unique<Patterns::Tuple>(
+          internal::default_map_separator[internal::RankInfo<T>::map_rank - 1],
           *Convert<Key>::to_pattern(),
-          *Convert<Value>::to_pattern(),
-          1,
-          1,
-          // We keep the same list separator of the previous level, as this is
-          // a map with only 1 possible entry
-          internal::default_list_separator[internal::RankInfo<T>::list_rank],
-          internal::default_map_separator[internal::RankInfo<T>::map_rank - 1]);
+          *Convert<Value>::to_pattern());
       }
 
       static std::string
@@ -2247,9 +2242,8 @@ namespace Patterns
                 const std::unique_ptr<Patterns::PatternBase> &pattern =
                   Convert<T>::to_pattern())
       {
-        std::unordered_map<Key, Value> m;
-        m.insert(t);
-        std::string s = Convert<decltype(m)>::to_string(m, pattern);
+        std::tuple<Key, Value> m(t);
+        std::string            s = Convert<decltype(m)>::to_string(m, pattern);
         AssertThrow(pattern->match(s), ExcNoMatch(s, pattern->description()));
         return s;
       }
@@ -2259,9 +2253,9 @@ namespace Patterns
                const std::unique_ptr<Patterns::PatternBase> &pattern =
                  Convert<T>::to_pattern())
       {
-        std::unordered_map<Key, Value> m;
+        std::tuple<Key, Value> m;
         m = Convert<decltype(m)>::to_value(s, pattern);
-        return *m.begin();
+        return m;
       }
     };
 
index b74fb5d1b54d69e0d868ee27e9417c85d5230772..339c7d83433d59e46298d13d87a1c177d53fc46f 100644 (file)
@@ -163,6 +163,150 @@ template <typename IndexType = boost::geometry::index::linear<16>,
 RTree<typename ContainerType::value_type, IndexType>
 pack_rtree(const ContainerType &container);
 
+/**
+ * Helper structure that allows one to extract a level from an RTree as a vector
+ * of BoundingBox objects.
+ *
+ * This structure implements a boost::geometry::index::detail::rtree::visitor
+ * object, allowing one to visit any existing RTree object, and return the
+ * vector of bounding boxes associated to a specific target level of the tree.
+ *
+ * Although possible, direct usage of this structure is cumbersome. The
+ * suggested usage of this class is through the helper function
+ * extract_rtree_level().
+ *
+ * @author Luca Heltai, 2020.
+ */
+template <typename Value,
+          typename Options,
+          typename Translator,
+          typename Box,
+          typename Allocators>
+struct ExtractLevelVisitor
+  : public boost::geometry::index::detail::rtree::visitor<
+      Value,
+      typename Options::parameters_type,
+      Box,
+      Allocators,
+      typename Options::node_tag,
+      true>::type
+{
+  /**
+   * Construct a vector @p boxes of BoundingBox objects corresponding to the
+   * @p target_level of the tree.
+   */
+  inline ExtractLevelVisitor(
+    Translator const & translator,
+    const unsigned int target_level,
+    std::vector<BoundingBox<boost::geometry::dimension<Box>::value>> &boxes);
+
+  /**
+   * An alias that identifies an InternalNode of the tree.
+   */
+  using InternalNode =
+    typename boost::geometry::index::detail::rtree::internal_node<
+      Value,
+      typename Options::parameters_type,
+      Box,
+      Allocators,
+      typename Options::node_tag>::type;
+
+  /**
+   * An alias that identifies a Leaf of the tree.
+   */
+  using Leaf = typename boost::geometry::index::detail::rtree::leaf<
+    Value,
+    typename Options::parameters_type,
+    Box,
+    Allocators,
+    typename Options::node_tag>::type;
+
+  /**
+   * Implements the visitor interface for InternalNode objects. If the node
+   * belongs to the @p target_leve, then fill the bounding box vector.
+   */
+  inline void
+  operator()(InternalNode const &node);
+
+  /**
+   * Implements the visitor interface for Leaf objects.
+   */
+  inline void
+  operator()(Leaf const &);
+
+  /**
+   * Translator interface, required by the boost implementation of the rtree.
+   */
+  Translator const &translator;
+
+  /**
+   * Store the level we are currently visiting.
+   */
+  size_t level;
+
+  /**
+   * The level we want to extract from the RTree object.
+   */
+  const size_t target_level;
+
+  /**
+   * A reference to the input vector of BoundingBox objects.
+   */
+  std::vector<BoundingBox<boost::geometry::dimension<Box>::value>> &boxes;
+};
+
+/**
+ * Given a RTree object @p rtree, and a target level @p level, return a vector
+ * of BoundingBox objects containing all the bounding boxes that make the given
+ * @p level of the @p rtree. This function is a convenient wrapper around the
+ * ExtractLevelVisitor class.
+ *
+ * Since an RTree object is a balanced tree, you can expect each entry of the
+ * resulting vector to contain roughly the same number of children, and
+ * ultimately, the same number of leaf objects. If you request for a level that
+ * is not present in the RTree, an empty vector is returned.
+ *
+ * A typical usage of this function is in the context of
+ * parallel::distributed::Triangulation objects, where one would like to
+ * construct a rough representation of the area which is covered by the locally
+ * owned cells of the active process, and exchange this information with other
+ * processes. The finest level of information is given by the leaves, which in
+ * this context would be the collection of all the bounding boxes associated
+ * to the locally owned cells of the triangulation. Exchanging this information
+ * with all participating processess would defeat the purpuse of parallel
+ * computations. If however one constructs an RTree containing these bounding
+ * boxes (for example, by calling
+ * GridTools::Cache::get_cell_bounding_boxes_rtree()), and then extracts one of
+ * the first levels of the RTree, only a handful of BoundingBoxes would be
+ * returned, allowing the user to have a very efficient description of the
+ * geometry of the domain, and of its distribution among processes.
+ *
+ * @author Luca Heltai, 2020.
+ */
+template <typename Rtree>
+inline std::vector<
+  BoundingBox<boost::geometry::dimension<typename Rtree::value_type>::value>>
+extract_rtree_level(Rtree const &tree, const unsigned int level)
+{
+  constexpr unsigned int dim =
+    boost::geometry::dimension<typename Rtree::value_type>::value;
+
+  using RtreeView =
+    boost::geometry::index::detail::rtree::utilities::view<Rtree>;
+  RtreeView rtv(tree);
+
+  std::vector<BoundingBox<dim>> boxes;
+
+  ExtractLevelVisitor<typename RtreeView::value_type,
+                      typename RtreeView::options_type,
+                      typename RtreeView::translator_type,
+                      typename RtreeView::box_type,
+                      typename RtreeView::allocators_type>
+    extract_level_visitor(rtv.translator(), level, boxes);
+  rtv.apply_visitor(extract_level_visitor);
+  return boxes;
+}
+
 
 
 // Inline and template functions
@@ -183,6 +327,80 @@ pack_rtree(const ContainerType &container)
   return pack_rtree<IndexType>(container.begin(), container.end());
 }
 
+
+
+template <typename Value,
+          typename Options,
+          typename Translator,
+          typename Box,
+          typename Allocators>
+ExtractLevelVisitor<Value, Options, Translator, Box, Allocators>::
+  ExtractLevelVisitor(
+    const Translator & translator,
+    const unsigned int target_level,
+    std::vector<BoundingBox<boost::geometry::dimension<Box>::value>> &boxes)
+  : translator(translator)
+  , level(0)
+  , target_level(target_level)
+  , boxes(boxes)
+{}
+
+
+
+template <typename Value,
+          typename Options,
+          typename Translator,
+          typename Box,
+          typename Allocators>
+void
+ExtractLevelVisitor<Value, Options, Translator, Box, Allocators>::
+operator()(const ExtractLevelVisitor::InternalNode &node)
+{
+  using ElmentsType =
+    typename boost::geometry::index::detail::rtree::elements_type<
+      InternalNode>::type;
+
+  const auto &elements = boost::geometry::index::detail::rtree::elements(node);
+
+  if (level == target_level)
+    {
+      const auto offset = boxes.size();
+      boxes.resize(offset + elements.size());
+
+      unsigned int i = offset;
+      for (typename ElmentsType::const_iterator it = elements.begin();
+           it != elements.end();
+           ++it)
+        {
+          boost::geometry::convert(it->first, boxes[i]);
+          ++i;
+        }
+      return;
+    }
+
+  const size_t level_backup = level;
+  ++level;
+
+  for (typename ElmentsType::const_iterator it = elements.begin();
+       it != elements.end();
+       ++it)
+    {
+      boost::geometry::index::detail::rtree::apply_visitor(*this, *it->second);
+    }
+
+  level = level_backup;
+}
+
+template <typename Value,
+          typename Options,
+          typename Translator,
+          typename Box,
+          typename Allocators>
+void
+ExtractLevelVisitor<Value, Options, Translator, Box, Allocators>::
+operator()(const ExtractLevelVisitor::Leaf &)
+{}
+
 #endif
 
 DEAL_II_NAMESPACE_CLOSE
diff --git a/tests/boost/extract_reprsentative_set_01.cc b/tests/boost/extract_reprsentative_set_01.cc
new file mode 100644 (file)
index 0000000..bab26d2
--- /dev/null
@@ -0,0 +1,49 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 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.
+//
+// ---------------------------------------------------------------------
+
+// Extract a representative vector of BoundingBox objects from an RTree
+
+#include <deal.II/base/patterns.h>
+
+#include <deal.II/boost_adaptors/bounding_box.h>
+#include <deal.II/boost_adaptors/point.h>
+
+#include <deal.II/numerics/rtree.h>
+
+#include <algorithm>
+
+#include "../tests.h"
+
+int
+main()
+{
+  initlog(1);
+
+  const unsigned int    N = 20;
+  std::vector<Point<2>> points(N);
+
+  for (auto &p : points)
+    p = random_point<2>();
+
+  // Make sure we have at most two points in each box
+  const auto tree = pack_rtree<boost::geometry::index::linear<2>>(points);
+
+  const auto boxes = extract_rtree_level(tree, 1);
+  deallog << "LEVEL 1:  N boxes: " << boxes.size() << std::endl;
+
+  for (const auto &b : boxes)
+    deallog << "Box: " << Patterns::Tools::to_string(b.get_boundary_points())
+            << std::endl;
+}
\ No newline at end of file
diff --git a/tests/boost/extract_reprsentative_set_01.output b/tests/boost/extract_reprsentative_set_01.output
new file mode 100644 (file)
index 0000000..6aa4bc3
--- /dev/null
@@ -0,0 +1,5 @@
+
+DEAL::LEVEL 1:  N boxes: 3
+DEAL::Box: 0.0163006, 0.108809 : 0.477397, 0.76823
+DEAL::Box: 0.137232, 0.771358 : 0.512932, 0.972775
+DEAL::Box: 0.61264, 0.197551 : 0.998925, 0.916195

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.