--- /dev/null
+New: GridTools::find_active_cell_around_point now allows you to specify an (optional) rtree, constructed
+from the used vertices of the triangulation. Once you have built a tree, querying for a nearest
+vertex is an O(log(N)) operation, where N is the number of used vertices. You can ask a GridTools::Cache
+object to return a tree that is compatible with the new function signature.
+The previous version of this function had a cost that was O(N^2) when the point was not in the cell_hint
+object.
+<br>
+(Luca Heltai, 2018/12/05)
+
# include <deal.II/lac/sparsity_tools.h>
+# include <deal.II/numerics/rtree.h>
+
# include <boost/archive/binary_iarchive.hpp>
# include <boost/archive/binary_oarchive.hpp>
# include <boost/optional.hpp>
* A version of the previous function that exploits an already existing
* map between vertices and cells, constructed using the function
* GridTools::vertex_to_cell_map, a map of vertex_to_cell_centers, obtained
- * through GridTools::vertex_to_cell_centers_directions, and a guess
- * `cell_hint`.
+ * through GridTools::vertex_to_cell_centers_directions, a guess
+ * `cell_hint`, and optionally an RTree constructed from the used
+ * vertices of the Triangulation. All of these structures can be queried
+ * from a GridTools::Cache object.
*
* @author Luca Heltai, Rene Gassmoeller, 2017
*/
const std::vector<std::vector<Tensor<1, spacedim>>> &vertex_to_cell_centers,
const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint =
typename MeshType<dim, spacedim>::active_cell_iterator(),
- const std::vector<bool> &marked_vertices = {});
+ const std::vector<bool> & marked_vertices = {},
+ const RTree<std::pair<Point<spacedim>, unsigned int>> &used_vertices_rtree =
+ {});
/**
* A version of the previous function where we use that mapping on a given
& vertex_to_cells,
const std::vector<std::vector<Tensor<1, spacedim>>> &vertex_to_cell_centers,
const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint,
- const std::vector<bool> &marked_vertices)
+ const std::vector<bool> & marked_vertices,
+ const RTree<std::pair<Point<spacedim>, unsigned int>> &used_vertices_rtree)
{
std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
}
else
{
- closest_vertex_index =
- GridTools::find_closest_vertex(mesh, p, marked_vertices);
+ if (!used_vertices_rtree.empty())
+ {
+ // If we have an rtree at our disposal, use it.
+ using ValueType = std::pair<Point<spacedim>, unsigned int>;
+ std::function<bool(const ValueType &)> marked;
+ if (marked_vertices.size() == mesh.n_vertices())
+ marked = [&marked_vertices](const ValueType &value) -> bool {
+ return marked_vertices[value.second];
+ };
+ else
+ marked = [](const ValueType &) -> bool { return true; };
+
+ std::vector<std::pair<Point<spacedim>, unsigned int>> res;
+ used_vertices_rtree.query(
+ boost::geometry::index::nearest(p, 1) &&
+ boost::geometry::index::satisfies(marked),
+ std::back_inserter(res));
+
+ // We should have one and only one result
+ AssertDimension(res.size(), 1);
+ closest_vertex_index = res[0].second;
+ }
+ else
+ {
+ closest_vertex_index =
+ GridTools::find_closest_vertex(mesh, p, marked_vertices);
+ }
vertex_to_point = p - mesh.get_vertices()[closest_vertex_index];
}
const auto &vertex_to_cells = cache.get_vertex_to_cell_map();
const auto &vertex_to_cell_centers =
cache.get_vertex_to_cell_centers_directions();
+ const auto &used_vertices_rtree = cache.get_used_vertices_rtree();
return find_active_cell_around_point(mapping,
mesh,
vertex_to_cells,
vertex_to_cell_centers,
cell_hint,
- marked_vertices);
+ marked_vertices,
+ used_vertices_rtree);
}
template <int spacedim>
const dealii::internal::ActiveCellIterator<deal_II_dimension,
deal_II_space_dimension,
X>::type &,
- const std::vector<bool> &);
+ const std::vector<bool> &,
+ const RTree<std::pair<Point<deal_II_space_dimension>, unsigned int>> &);
template std::vector<BoundingBox<deal_II_space_dimension>>
compute_mesh_predicate_bounding_box<X>(