--- /dev/null
+Changed: GridTools::distributed_compute_point_locations now takes as input the
+global description of the manifold using bounding boxes.
+<br>
+(Giovanni Alzetta, 2018/04/16)
+
* @param[in] local_points the array of points owned by the current process. Every
* process can have a different array of points which can be empty and not
* contained within the locally owned part of the triangulation
- * @param[in] local_bbox the description of the locally owned part of the mesh made
- * with bounding boxes. It can be obtained from
- * GridTools::compute_mesh_predicate_bounding_box
+ * @param[in] global_bboxes a vector of vectors of bounding boxes; it describes
+ * the locally owned part of the mesh for each process.
+ * The bounding boxes describing which part of the mesh is locally owned by
+ * process with rank rk are contained in global_bboxes[rk].
+ * The local description can be obtained from
+ * GridTools::compute_mesh_predicate_bounding_box ; then the global one can
+ * be obtained using either GridTools::exchange_local_bounding_boxes
+ * or Utilities::MPI::all_gather
* @return A tuple containing the quadrature information
*
* The elements of the output tuple are:
return_type
#endif
distributed_compute_point_locations
- (const GridTools::Cache<dim,spacedim> &cache,
- const std::vector<Point<spacedim> > &local_points,
- const std::vector< BoundingBox<spacedim> > &local_bbox);
+ (const GridTools::Cache<dim,spacedim> &cache,
+ const std::vector<Point<spacedim> > &local_points,
+ const std::vector< std::vector< BoundingBox<spacedim> > > &global_bboxes);
/**
* Return a map of index:Point<spacedim>, containing the used vertices of the
std::vector< std::vector< unsigned int > >
>
distributed_compute_point_locations
- (const GridTools::Cache<dim,spacedim> &cache,
- const std::vector<Point<spacedim> > &local_points,
- const std::vector< BoundingBox<spacedim> > &local_bbox)
+ (const GridTools::Cache<dim,spacedim> &cache,
+ const std::vector<Point<spacedim> > &local_points,
+ const std::vector< std::vector< BoundingBox<spacedim> > > &global_bboxes)
{
#ifndef DEAL_II_WITH_MPI
(void)cache;
(void)local_points;
- (void)local_bbox;
+ (void)global_bboxes;
Assert(false, ExcMessage("GridTools::distributed_compute_point_locations() requires MPI."));
std::tuple<
std::vector< typename Triangulation<dim, spacedim>::active_cell_iterator >,
internal::distributed_cptloc::cell_hash<dim,spacedim> >
temporary_unmap;
- // Obtaining the global mesh description through an all to all communication
- std::vector< std::vector< BoundingBox<spacedim> > > global_bounding_boxes;
- global_bounding_boxes = Utilities::MPI::all_gather(mpi_communicator,local_bbox);
-
// Step 1 (part 1): Using the bounding boxes to guess the owner of each points
// in local_points
unsigned int my_rank = Utilities::MPI::this_mpi_process(mpi_communicator);
std::tuple< std::vector< std::vector< unsigned int > >, std::map< unsigned int, unsigned int >,
std::map< unsigned int, std::vector< unsigned int > > > guessed_points;
guessed_points =
- GridTools::guess_point_owner(global_bounding_boxes, local_points);
+ GridTools::guess_point_owner(global_bboxes, local_points);
// Preparing to call compute point locations on points which are/might be
// local
std::vector< std::vector< unsigned int > >
>
distributed_compute_point_locations
- (const Cache< deal_II_dimension, deal_II_space_dimension > &,
- const std::vector< Point< deal_II_space_dimension > > &,
- const std::vector< BoundingBox< deal_II_space_dimension > > &);
+ (const Cache< deal_II_dimension, deal_II_space_dimension > &,
+ const std::vector< Point< deal_II_space_dimension > > &,
+ const std::vector< std::vector< BoundingBox< deal_II_space_dimension > > > &);
\}
#endif
std::vector< BoundingBox<dim> > local_bbox = GridTools::compute_mesh_predicate_bounding_box
(cache.get_triangulation(),
std::function<bool (const typename Triangulation<dim>::active_cell_iterator &)>(locally_owned_cell_predicate),
- 1, true, 4); // These options should be passed
- // Using the distributed version of compute point location
+ 1, true, 4);
+
+ // Obtaining the global mesh description through an all to all communication
+ std::vector< std::vector< BoundingBox<dim> > > global_bboxes;
+ global_bboxes = Utilities::MPI::all_gather(mpi_communicator,local_bbox);
- // Using the distributed version
- auto output_tuple = distributed_compute_point_locations(cache,points,local_bbox);
+ // Using the distributed version of compute point location
+ auto output_tuple = distributed_compute_point_locations(cache,points,global_bboxes);
// Testing in serial against the serial version
auto cell_qpoint_map = GridTools::compute_point_locations(cache,points);
(cache.get_triangulation(), locally_owned_cell_predicate,
1, true, 4);
+ // Obtaining the global mesh description through an all to all communication
+ std::vector< std::vector< BoundingBox<dim> > > global_bboxes;
+ global_bboxes = Utilities::MPI::all_gather(mpi_communicator,local_bbox);
+
// Using the distributed version of compute point location
auto output_tuple = distributed_compute_point_locations
- (cache,loc_owned_points,local_bbox);
+ (cache,loc_owned_points,global_bboxes);
deallog << "Comparing results" << std::endl;
const auto &output_cells = std::get<0>(output_tuple);
const auto &output_qpoints = std::get<1>(output_tuple);