--- /dev/null
+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)
+
# 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>
# 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>
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
#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 */
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);
}
#if deal_II_dimension <= deal_II_space_dimension
namespace GridTools
\{
-
template double
diameter(
const Triangulation<deal_II_dimension, deal_II_space_dimension> &);
&distorted_cells,
Triangulation<deal_II_dimension, deal_II_space_dimension>
&triangulation);
-
# endif
-
\}
#endif
}